Generate and return a resized and centered copy of the input array
within the newly allocated arrayCentered
and filling the newly added elements (if any) with the user-specified fill
.
Additionally, if left and right margins lmsize, rmsize
are specified, resize the input allocatable array(lb:ub)
to arrayCentered(lb : lb + size + lmsize + rmsize - 1)
while centering the original array
contents within the range (lb + lmsize : lb + lmsize + size)
in the newly allocated array.
- Parameters
-
[in,out] | array | : The input allocatable array of shape (:) of either
-
type
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU), or
-
type
integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64), or
-
type
logical of kind any supported by the processor (e.g., LK), 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).
On output, the array will be reallocated to the requested new size, with the same lower bound as before. |
[in] | size | : The input scalar integer of default kind IK representing the length of the section of the array within which the array contents will be centered.
If the lmsize and rmsize input arguments are missing, then size represents the total length of the output array, otherwise the length of the output array is lmsize + size + rmsize . |
[in] | lmsize | : The input scalar integer of default kind IK representing the size of the left-margin of the output array
(optional, default = 0 . It can be present only if the input argument rmsize is also present.) |
[in] | lmsize | : The input scalar integer of default kind IK representing the size of the right-margin of the output array
(optional, default = 0 . It can be present only if the input argument lmsize is also present.) |
[in] | fill | : The input scalar of the same type and kind as the input array containing the value to fill the new elements (if any) of array surrounding the original array contents (excluding the margins). If array is of type character , then
-
the equality
len(fill) == 1 must also hold if array is a scalar string.
-
the equality
len(fill) == len(array) must also hold if array is an array of strings.
(optional, if missing, the new elements will remain uninitialized). |
[in] | lmfill | : The input scalar of the same type and kind as the input array containing the value to fill the left margin (if any) of newly allocated array . If array is of type character , then
-
the equality
len(lmfill) == 1 must also hold if array is a scalar string.
-
the equality
len(lmfill) == len(array) must also hold if array is an array of strings.
(optional, if missing, the left-margin will remain uninitialized. It can be present only if lmsize and rmsize are also present.) |
[in] | rmfill | : The input scalar of the same type and kind as the input array containing the value to fill the right margin (if any) of newly allocated array . If array is of type character , then
-
the equality
len(rmfill) == 1 must also hold if array is a scalar string.
-
the equality
len(rmfill) == len(array) must also hold if array is an array of strings.
(optional, if missing, the right-margin will remain uninitialized. It can be present only if lmsize and rmsize are also present.) |
- Returns
arrayCentered
: The output object of the same type, kind, and rank as the input array
whose size is size + lmsize + rmsize
, whose contents are the same as the contents of array
but centered in the range (1 + lmsize : lmsize + size)
while any new elements are optionally filled with fill
and the sections (1 : lmsize)
and (lmsize + size + 1 : lmsize + size + rmsize)
corresponding to the array margins are filled with lmfill
and rmfill
respectively.
Possible calling interfaces ⛓
arrayCentered
= getCentered(array, size, lmsize, rmsize, fill
= fill, lmfill
= lmfill, rmfill
= rmfill)
!
Generate and return a resized and centered copy of the input array within the newly allocated arrayCe...
This module contains procedures and generic interfaces for resizing an input array and centering the ...
- Warning
- Note that the new elements of the newly allocated output
arrayCentered
are not initialized to any particular value if fill
is missing.
Therefore, the contents of the new elements of the output arrayCentered
is processor dependent, frequently meaningless and should not be relied upon.
However, if fill
is specified, then all new elements will be filled with the specified fill
value.
-
Similarly, the margin elements will not be initialized to any particular values unless the corresponding
lmfill
or rmfill
arguments are present.
-
When the specified new array size is smaller than the original, the corresponding elements of the old array will be also trimmed from the resized shrunk array.
In such a case, there will be no new elements to initialize their values and fill
will have no effects on the output.
-
The input
size
, lmsize
, rmsize
arguments must be non-negative input arguments.
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
- The sole purpose of this generic interface is to provide a convenient but fast method of generating a resized and centered copy of an array.
See pm_arrayResize and pm_arrayCenter for the relevant benchmarks.
- See also
- setResized
setRebound
setRefilled
setRebilled
getCoreHalo
setCoreHalo
getCentered
setCentered
getPadded
setPadded
isFailedGetShellWidth
isFailedGetShellHeight
Example usage ⛓
13 character(:, SK),
allocatable :: string_SK
14 character(
2, SK),
allocatable :: array_SK(:)
15 logical(LK) ,
allocatable :: array_LK(:)
16 integer(IK) ,
allocatable :: array_IK(:)
17 complex(CK) ,
allocatable :: array_CK(:)
18 real(RK) ,
allocatable :: array_RK(:)
20 type(display_type) :: disp
25 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show(
"! Expand and center an array to a new larger size.")
28 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
34 allocate(array_SK(
3:
8), source
= [SK_
"AA", SK_
"BB", SK_
"CC", SK_
"DD", SK_
"EE", SK_
"FF"])
35 allocate(array_IK(
3:
8),
source = [
1_IK,
2_IK,
3_IK,
4_IK,
5_IK,
6_IK])
36 allocate(array_RK(
3:
8),
source = [
1._RK,
2._RK,
3._RK,
4._RK,
5._RK,
6._RK])
37 allocate(array_CK(
3:
8),
source = [(
1._CK,
-1._CK), (
2._CK,
-2._CK), (
3._CK,
-3._CK), (
4._CK,
-4._CK), (
5._CK,
-5._CK), (
6._CK,
-6._CK)])
38 allocate(array_LK(
3:
8),
source = [
.true._LK,
.true._LK,
.true._LK,
.true._LK,
.true._LK,
.true._LK])
41 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
42 call disp%show(
"! Center character scalar.")
43 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
46 string_SK
= SK_
"ABCDEF"
48 call disp%show( string_SK, deliml
= SK_
"""" )
51 call disp%show(
"string_SK = getCentered(string_SK, size = 8_IK, fill = SK_'-')")
52 string_SK
= getCentered(string_SK, size
= 8_IK, fill
= SK_
'-')
54 call disp%show( string_SK, deliml
= SK_
"""" )
59 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
60 call disp%show(
"! Center character array.")
61 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
65 call disp%show( array_SK, deliml
= SK_
"""" )
66 call disp%show(
"lbound(array_SK), ubound(array_SK)")
67 call disp%show([
lbound(array_SK),
ubound(array_SK)])
68 call disp%show(
"array_SK = getCentered(array_SK, size = 8_IK, fill = SK_'%%')")
69 array_SK
= getCentered(array_SK, size
= 8_IK, fill
= SK_
'%%')
71 call disp%show( array_SK, deliml
= SK_
"""" )
72 call disp%show(
"lbound(array_SK), ubound(array_SK)")
73 call disp%show([
lbound(array_SK),
ubound(array_SK)])
76 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
77 call disp%show(
"! Center logical array.")
78 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
83 call disp%show(
"lbound(array_LK), ubound(array_LK)")
84 call disp%show([
lbound(array_LK),
ubound(array_LK)])
85 call disp%show(
"array_LK = getCentered(array_LK, size = 8_IK, fill = .false._LK)")
86 array_LK
= getCentered(array_LK,
size = 8_IK, fill
= .false._LK)
89 call disp%show(
"lbound(array_LK), ubound(array_LK)")
90 call disp%show([
lbound(array_LK),
ubound(array_LK)])
93 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
94 call disp%show(
"! Center integer array.")
95 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
100 call disp%show(
"lbound(array_IK), ubound(array_IK)")
101 call disp%show([
lbound(array_IK),
ubound(array_IK)])
102 call disp%show(
"array_IK = getCentered(array_IK, size = 8_IK, fill = -9999_IK)")
103 array_IK
= getCentered(array_IK,
size = 8_IK, fill
= -9999_IK)
106 call disp%show(
"lbound(array_IK), ubound(array_IK)")
107 call disp%show([
lbound(array_IK),
ubound(array_IK)])
110 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
111 call disp%show(
"! Center complex array.")
112 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
117 call disp%show(
"lbound(array_CK), ubound(array_CK)")
118 call disp%show([
lbound(array_CK),
ubound(array_CK)])
119 call disp%show(
"array_CK = getCentered(array_CK, size = 8_IK, fill = (-999._CK,-999._CK))")
120 array_CK
= getCentered(array_CK,
size = 8_IK, fill
= (
-999._CK,
-999._CK))
123 call disp%show(
"lbound(array_CK), ubound(array_CK)")
124 call disp%show([
lbound(array_CK),
ubound(array_CK)])
127 call disp%show(
"!%%%%%%%%%%%%%%%%%%%")
128 call disp%show(
"! Center real array.")
129 call disp%show(
"!%%%%%%%%%%%%%%%%%%%")
134 call disp%show(
"lbound(array_RK), ubound(array_RK)")
135 call disp%show([
lbound(array_RK),
ubound(array_RK)])
136 call disp%show(
"array_RK = getCentered(array_RK, size = 8_IK, fill = 0._RK)")
137 array_RK
= getCentered(array_RK,
size = 8_IK, fill
= 0._RK)
140 call disp%show(
"lbound(array_RK), ubound(array_RK)")
141 call disp%show([
lbound(array_RK),
ubound(array_RK)])
145 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
146 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
147 call disp%show(
"! Shrink an array to a new smaller size. Note that fill has no effect because the array shrinks.")
148 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
149 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
154 string_SK
= SK_
"ABCDEF"
155 allocate(array_SK(
3:
8), source
= [SK_
"AA", SK_
"BB", SK_
"CC", SK_
"DD", SK_
"EE", SK_
"FF"])
156 allocate(array_IK(
3:
8),
source = [
1_IK,
2_IK,
3_IK,
4_IK,
5_IK,
6_IK])
157 allocate(array_RK(
3:
8),
source = [
1._RK,
2._RK,
3._RK,
4._RK,
5._RK,
6._RK])
158 allocate(array_CK(
3:
8),
source = [(
1._CK,
-1._CK), (
2._CK,
-2._CK), (
3._CK,
-3._CK), (
4._CK,
-4._CK), (
5._CK,
-5._CK), (
6._CK,
-6._CK)])
159 allocate(array_LK(
3:
8),
source = [
.true._LK,
.true._LK,
.true._LK,
.true._LK,
.true._LK,
.true._LK])
162 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
163 call disp%show(
"! Center character scalar.")
164 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
168 call disp%show( string_SK, deliml
= SK_
"""" )
171 call disp%show(
"string_SK = getCentered(string_SK, size = 3_IK, fill = '~')")
172 string_SK
= getCentered(string_SK, size
= 3_IK, fill
= '~')
174 call disp%show( string_SK, deliml
= SK_
"""" )
179 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
180 call disp%show(
"! Center character array.")
181 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
185 call disp%show( array_SK, deliml
= SK_
"""" )
186 call disp%show(
"lbound(array_SK), ubound(array_SK)")
187 call disp%show([
lbound(array_SK),
ubound(array_SK)])
188 call disp%show(
"array_SK = getCentered(array_SK, size = 3_IK, fill = '**')")
189 array_SK
= getCentered(array_SK, size
= 3_IK, fill
= '**')
191 call disp%show( array_SK, deliml
= SK_
"""" )
192 call disp%show(
"lbound(array_SK), ubound(array_SK)")
193 call disp%show([
lbound(array_SK),
ubound(array_SK)])
196 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
197 call disp%show(
"! Center logical array.")
198 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
203 call disp%show(
"lbound(array_LK), ubound(array_LK)")
204 call disp%show([
lbound(array_LK),
ubound(array_LK)])
205 call disp%show(
"array_LK = getCentered(array_LK, size = 3_IK, fill = .false._LK)")
206 array_LK
= getCentered(array_LK,
size = 3_IK, fill
= .false._LK)
209 call disp%show(
"lbound(array_LK), ubound(array_LK)")
210 call disp%show([
lbound(array_LK),
ubound(array_LK)])
213 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
214 call disp%show(
"! Center integer array.")
215 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
220 call disp%show(
"lbound(array_IK), ubound(array_IK)")
221 call disp%show([
lbound(array_IK),
ubound(array_IK)])
222 call disp%show(
"array_IK = getCentered(array_IK, size = 3_IK, fill = 0_IK)")
223 array_IK
= getCentered(array_IK,
size = 3_IK, fill
= 0_IK)
226 call disp%show(
"lbound(array_IK), ubound(array_IK)")
227 call disp%show([
lbound(array_IK),
ubound(array_IK)])
230 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
231 call disp%show(
"! Center complex array.")
232 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
237 call disp%show(
"lbound(array_CK), ubound(array_CK)")
238 call disp%show([
lbound(array_CK),
ubound(array_CK)])
239 call disp%show(
"array_CK = getCentered(array_CK, size = 3_IK, fill = (0._CK, 0._CK))")
240 array_CK
= getCentered(array_CK,
size = 3_IK, fill
= (
0._CK,
0._CK))
243 call disp%show(
"lbound(array_CK), ubound(array_CK)")
244 call disp%show([
lbound(array_CK),
ubound(array_CK)])
247 call disp%show(
"!%%%%%%%%%%%%%%%%%%%")
248 call disp%show(
"! Center real array.")
249 call disp%show(
"!%%%%%%%%%%%%%%%%%%%")
254 call disp%show(
"lbound(array_RK), ubound(array_RK)")
255 call disp%show([
lbound(array_RK),
ubound(array_RK)])
256 call disp%show(
"array_RK = getCentered(array_RK, size = 3_IK, fill = +0._RK)")
257 array_RK
= getCentered(array_RK,
size = 3_IK, fill
= +0._RK)
260 call disp%show(
"lbound(array_RK), ubound(array_RK)")
261 call disp%show([
lbound(array_RK),
ubound(array_RK)])
264 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
265 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
266 call disp%show(
"! Example of perfect array centering.")
267 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
268 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
273 call disp%show( string_SK, deliml
= SK_
"""" )
276 call disp%show(
"string_SK = getCentered(string_SK, size = 7_IK, fill = SK_'-')")
277 string_SK
= getCentered(string_SK, size
= 7_IK, fill
= SK_
'-')
279 call disp%show( string_SK, deliml
= SK_
"""" )
286 call disp%show( string_SK, deliml
= SK_
"""" )
289 call disp%show(
"string_SK = getCentered(string_SK, size = 1_IK, fill = SK_'-')")
290 string_SK
= getCentered(string_SK, size
= 1_IK, fill
= SK_
'-')
292 call disp%show( string_SK, deliml
= SK_
"""" )
299 call disp%show( string_SK, deliml
= SK_
"""" )
302 call disp%show(
"string_SK = getCentered(string_SK, size = 6_IK, fill = SK_'-')")
303 string_SK
= getCentered(string_SK, size
= 6_IK, fill
= SK_
'-')
305 call disp%show( string_SK, deliml
= SK_
"""" )
311 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
312 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
313 call disp%show(
"! Example of imperfect array centering.")
314 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
315 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
320 call disp%show( string_SK, deliml
= SK_
"""" )
323 call disp%show(
"string_SK = getCentered(string_SK, size = 6_IK, fill = SK_'-')")
324 string_SK
= getCentered(string_SK, size
= 6_IK, fill
= SK_
'-')
326 call disp%show( string_SK, deliml
= SK_
"""" )
333 call disp%show( string_SK, deliml
= SK_
"""" )
336 call disp%show(
"string_SK = getCentered(string_SK, size = 2_IK, fill = SK_'-')")
337 string_SK
= getCentered(string_SK, size
= 2_IK, fill
= SK_
'-')
339 call disp%show( string_SK, deliml
= SK_
"""" )
346 call disp%show( string_SK, deliml
= SK_
"""" )
349 call disp%show(
"string_SK = getCentered(string_SK, size = 3_IK, fill = SK_'-')")
350 string_SK
= getCentered(string_SK, size
= 3_IK, fill
= SK_
'-')
352 call disp%show( string_SK, deliml
= SK_
"""" )
359 call disp%show( string_SK, deliml
= SK_
"""" )
362 call disp%show(
"string_SK = getCentered(string_SK, size = 1_IK, fill = SK_'-')")
363 string_SK
= getCentered(string_SK, size
= 1_IK, fill
= SK_
'-')
365 call disp%show( string_SK, deliml
= SK_
"""" )
371 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
372 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
373 call disp%show(
"! Centering array contents within a given size, padded by left and right margins.")
374 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
375 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
378 string_SK
= SK_
"In God We Trust"
380 call disp%show( string_SK , deliml
= SK_
"""" )
383 call disp%show(
"string_SK = getCentered(string_SK, size = 45_IK, fill = SK_' ', lmsize = 4_IK, rmsize = 4_IK, lmfill = SK_'*', rmfill = SK_'*')")
384 string_SK
= getCentered(string_SK, size
= 45_IK, fill
= SK_
' ', lmsize
= 4_IK, rmsize
= 4_IK, lmfill
= SK_
'*', rmfill
= SK_
'*')
386 call disp%show( string_SK, deliml
= SK_
"""" )
391 string_SK
= SK_
"ABCDEF"
393 call disp%show( string_SK, deliml
= SK_
"""" )
396 call disp%show(
"string_SK = getCentered(string_SK, size = 8_IK, lmsize = 3_IK, rmsize = 5_IK, fill = SK_'-', lmfill = SK_'*', rmfill = SK_'+')")
397 string_SK
= getCentered(string_SK, size
= 8_IK, lmsize
= 3_IK, rmsize
= 5_IK, fill
= SK_
'-', lmfill
= SK_
'*', rmfill
= SK_
'+')
399 call disp%show( string_SK, deliml
= SK_
"""" )
404 string_SK
= SK_
"ABCDEF"
406 call disp%show( string_SK, deliml
= SK_
"""" )
409 call disp%show(
"string_SK = getCentered(string_SK, size = 2_IK, lmsize = 3_IK, rmsize = 5_IK, fill = SK_'-', lmfill = SK_'*', rmfill = SK_'+')")
410 string_SK
= getCentered(string_SK, size
= 2_IK, lmsize
= 3_IK, rmsize
= 5_IK, fill
= SK_
'-', lmfill
= SK_
'*', rmfill
= SK_
'+')
412 call disp%show( string_SK, deliml
= SK_
"""" )
420 if (
allocated(array_SK))
deallocate(array_SK)
421 if (
allocated(array_IK))
deallocate(array_IK)
422 if (
allocated(array_RK))
deallocate(array_RK)
423 if (
allocated(array_CK))
deallocate(array_CK)
424 if (
allocated(array_LK))
deallocate(array_LK)
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 ⛓
17string_SK
= getCentered(string_SK, size
= 8_IK, fill
= SK_
'-')
28"AA",
"BB",
"CC",
"DD",
"EE",
"FF"
29lbound(array_SK),
ubound(array_SK)
31array_SK
= getCentered(array_SK, size
= 8_IK, fill
= SK_
'%%')
33"%%",
"AA",
"BB",
"CC",
"DD",
"EE",
"FF",
"%%"
34lbound(array_SK),
ubound(array_SK)
43lbound(array_LK),
ubound(array_LK)
45array_LK
= getCentered(array_LK,
size = 8_IK, fill
= .false._LK)
48lbound(array_LK),
ubound(array_LK)
57lbound(array_IK),
ubound(array_IK)
59array_IK
= getCentered(array_IK,
size = 8_IK, fill
= -9999_IK)
61-9999,
+1,
+2,
+3,
+4,
+5,
+6,
-9999
62lbound(array_IK),
ubound(array_IK)
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)
71lbound(array_CK),
ubound(array_CK)
73array_CK
= getCentered(array_CK,
size = 8_IK, fill
= (
-999._CK,
-999._CK))
75(
-999.00000000000000,
-999.00000000000000), (
+1.0000000000000000,
-1.0000000000000000), (
+2.0000000000000000,
-2.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000), (
+5.0000000000000000,
-5.0000000000000000), (
+6.0000000000000000,
-6.0000000000000000), (
-999.00000000000000,
-999.00000000000000)
76lbound(array_CK),
ubound(array_CK)
84+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
85lbound(array_RK),
ubound(array_RK)
87array_RK
= getCentered(array_RK,
size = 8_IK, fill
= 0._RK)
89+0.0000000000000000,
+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000,
+0.0000000000000000
90lbound(array_RK),
ubound(array_RK)
108string_SK
= getCentered(string_SK, size
= 3_IK, fill
= '~')
119"AA",
"BB",
"CC",
"DD",
"EE",
"FF"
120lbound(array_SK),
ubound(array_SK)
122array_SK
= getCentered(array_SK, size
= 3_IK, fill
= '**')
125lbound(array_SK),
ubound(array_SK)
134lbound(array_LK),
ubound(array_LK)
136array_LK
= getCentered(array_LK,
size = 3_IK, fill
= .false._LK)
139lbound(array_LK),
ubound(array_LK)
147+1,
+2,
+3,
+4,
+5,
+6
148lbound(array_IK),
ubound(array_IK)
150array_IK
= getCentered(array_IK,
size = 3_IK, fill
= 0_IK)
153lbound(array_IK),
ubound(array_IK)
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)
162lbound(array_CK),
ubound(array_CK)
164array_CK
= getCentered(array_CK,
size = 3_IK, fill
= (
0._CK,
0._CK))
166(
+2.0000000000000000,
-2.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000)
167lbound(array_CK),
ubound(array_CK)
175+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
176lbound(array_RK),
ubound(array_RK)
178array_RK
= getCentered(array_RK,
size = 3_IK, fill
= +0._RK)
180+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000
181lbound(array_RK),
ubound(array_RK)
194string_SK
= getCentered(string_SK, size
= 7_IK, fill
= SK_
'-')
204string_SK
= getCentered(string_SK, size
= 1_IK, fill
= SK_
'-')
214string_SK
= getCentered(string_SK, size
= 6_IK, fill
= SK_
'-')
231string_SK
= getCentered(string_SK, size
= 6_IK, fill
= SK_
'-')
241string_SK
= getCentered(string_SK, size
= 2_IK, fill
= SK_
'-')
251string_SK
= getCentered(string_SK, size
= 3_IK, fill
= SK_
'-')
261string_SK
= getCentered(string_SK, size
= 1_IK, fill
= SK_
'-')
278string_SK
= getCentered(string_SK, size
= 45_IK, fill
= SK_
' ', lmsize
= 4_IK, rmsize
= 4_IK, lmfill
= SK_
'*', rmfill
= SK_
'*')
280"**** In God We Trust ****"
288string_SK
= getCentered(string_SK, size
= 8_IK, lmsize
= 3_IK, rmsize
= 5_IK, fill
= SK_
'-', lmfill
= SK_
'*', rmfill
= SK_
'+')
298string_SK
= getCentered(string_SK, size
= 2_IK, lmsize
= 3_IK, rmsize
= 5_IK, fill
= SK_
'-', lmfill
= SK_
'*', rmfill
= SK_
'+')
- Test:
- test_pm_arrayCenter
- Todo:
- Normal Priority: Two new optional input scalar
lbcold
and ubcold
arguments can be added to procedures to specify a subset of the contents of the original array that has to be kept in the newly allocated centered array.
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 Austin
Definition at line 169 of file pm_arrayCenter.F90.