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

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.
More...

Detailed Description

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
  1. the equality len(fill) == 1 must also hold if array is a scalar string.
  2. 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
  1. the equality len(lmfill) == 1 must also hold if array is a scalar string.
  2. 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
  1. the equality len(rmfill) == 1 must also hold if array is a scalar string.
  2. 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, fill = fill) ! resize to `(lb:lb+size-1)`, center contents within the new limits, and optionally fill the new elements with `fill`.
arrayCentered = getCentered(array, size, lmsize, rmsize, fill = fill, lmfill = lmfill, rmfill = rmfill) ! resize to `(lb:lb+size+lmsize+rmsize-1)`, center contents, fill the new elements, left and right margins.
!
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

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 character(:, SK), allocatable :: string_SK
14 character(2, SK), allocatable :: array_SK(:) ! Can be any processor-supported kind.
15 logical(LK) , allocatable :: array_LK(:) ! Can be any processor-supported kind.
16 integer(IK) , allocatable :: array_IK(:) ! Can be any processor-supported kind.
17 complex(CK) , allocatable :: array_CK(:) ! Can be any processor-supported kind.
18 real(RK) , allocatable :: array_RK(:) ! Can be any processor-supported kind.
19
20 type(display_type) :: disp
21
22 disp = display_type(file = "main.out.F90")
23
24 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%skip()
31
32 call reset()
33
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])
39
40 call disp%skip()
41 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
42 call disp%show("! Center character scalar.")
43 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
44 call disp%skip()
45
46 string_SK = SK_"ABCDEF"
47 call disp%show("string_SK")
48 call disp%show( string_SK, deliml = SK_"""" )
49 call disp%show("len(string_SK)")
50 call disp%show( len(string_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_'-')
53 call disp%show("string_SK")
54 call disp%show( string_SK, deliml = SK_"""" )
55 call disp%show("len(string_SK)")
56 call disp%show( len(string_SK) )
57
58 call disp%skip()
59 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
60 call disp%show("! Center character array.")
61 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
62 call disp%skip()
63
64 call disp%show("array_SK")
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_'%%')
70 call disp%show("array_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)])
74
75 call disp%skip()
76 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
77 call disp%show("! Center logical array.")
78 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
79 call disp%skip()
80
81 call disp%show("array_LK")
82 call disp%show( array_LK )
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)
87 call disp%show("array_LK")
88 call disp%show( array_LK )
89 call disp%show("lbound(array_LK), ubound(array_LK)")
90 call disp%show([lbound(array_LK), ubound(array_LK)])
91
92 call disp%skip()
93 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
94 call disp%show("! Center integer array.")
95 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
96 call disp%skip()
97
98 call disp%show("array_IK")
99 call disp%show( array_IK )
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)
104 call disp%show("array_IK")
105 call disp%show( array_IK )
106 call disp%show("lbound(array_IK), ubound(array_IK)")
107 call disp%show([lbound(array_IK), ubound(array_IK)])
108
109 call disp%skip()
110 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
111 call disp%show("! Center complex array.")
112 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
113 call disp%skip()
114
115 call disp%show("array_CK")
116 call disp%show( array_CK )
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))
121 call disp%show("array_CK")
122 call disp%show( array_CK )
123 call disp%show("lbound(array_CK), ubound(array_CK)")
124 call disp%show([lbound(array_CK), ubound(array_CK)])
125
126 call disp%skip()
127 call disp%show("!%%%%%%%%%%%%%%%%%%%")
128 call disp%show("! Center real array.")
129 call disp%show("!%%%%%%%%%%%%%%%%%%%")
130 call disp%skip()
131
132 call disp%show("array_RK")
133 call disp%show( array_RK )
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)
138 call disp%show("array_RK")
139 call disp%show( array_RK )
140 call disp%show("lbound(array_RK), ubound(array_RK)")
141 call disp%show([lbound(array_RK), ubound(array_RK)])
142
143
144 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
150 call disp%skip()
151
152 call reset()
153
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])
160
161 call disp%skip()
162 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
163 call disp%show("! Center character scalar.")
164 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
165 call disp%skip()
166
167 call disp%show("string_SK")
168 call disp%show( string_SK, deliml = SK_"""" )
169 call disp%show("len(string_SK)")
170 call disp%show( len(string_SK) )
171 call disp%show("string_SK = getCentered(string_SK, size = 3_IK, fill = '~')")
172 string_SK = getCentered(string_SK, size = 3_IK, fill = '~')
173 call disp%show("string_SK")
174 call disp%show( string_SK, deliml = SK_"""" )
175 call disp%show("len(string_SK)")
176 call disp%show( len(string_SK) )
177
178 call disp%skip()
179 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
180 call disp%show("! Center character array.")
181 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
182 call disp%skip()
183
184 call disp%show("array_SK")
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 = '**')
190 call disp%show("array_SK")
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)])
194
195 call disp%skip()
196 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
197 call disp%show("! Center logical array.")
198 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
199 call disp%skip()
200
201 call disp%show("array_LK")
202 call disp%show( array_LK )
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)
207 call disp%show("array_LK")
208 call disp%show( array_LK )
209 call disp%show("lbound(array_LK), ubound(array_LK)")
210 call disp%show([lbound(array_LK), ubound(array_LK)])
211
212 call disp%skip()
213 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
214 call disp%show("! Center integer array.")
215 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
216 call disp%skip()
217
218 call disp%show("array_IK")
219 call disp%show( array_IK )
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)
224 call disp%show("array_IK")
225 call disp%show( array_IK )
226 call disp%show("lbound(array_IK), ubound(array_IK)")
227 call disp%show([lbound(array_IK), ubound(array_IK)])
228
229 call disp%skip()
230 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
231 call disp%show("! Center complex array.")
232 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
233 call disp%skip()
234
235 call disp%show("array_CK")
236 call disp%show( array_CK )
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))
241 call disp%show("array_CK")
242 call disp%show( array_CK )
243 call disp%show("lbound(array_CK), ubound(array_CK)")
244 call disp%show([lbound(array_CK), ubound(array_CK)])
245
246 call disp%skip()
247 call disp%show("!%%%%%%%%%%%%%%%%%%%")
248 call disp%show("! Center real array.")
249 call disp%show("!%%%%%%%%%%%%%%%%%%%")
250 call disp%skip()
251
252 call disp%show("array_RK")
253 call disp%show( array_RK )
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)
258 call disp%show("array_RK")
259 call disp%show( array_RK )
260 call disp%show("lbound(array_RK), ubound(array_RK)")
261 call disp%show([lbound(array_RK), ubound(array_RK)])
262
263 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
269 call disp%skip()
270
271 string_SK = SK_"ABC"
272 call disp%show("string_SK")
273 call disp%show( string_SK, deliml = SK_"""" )
274 call disp%show("len(string_SK)")
275 call disp%show( len(string_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_'-')
278 call disp%show("string_SK")
279 call disp%show( string_SK, deliml = SK_"""" )
280 call disp%show("len(string_SK)")
281 call disp%show( len(string_SK) )
282 call disp%skip()
283
284 string_SK = SK_"ABC"
285 call disp%show("string_SK")
286 call disp%show( string_SK, deliml = SK_"""" )
287 call disp%show("len(string_SK)")
288 call disp%show( len(string_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_'-')
291 call disp%show("string_SK")
292 call disp%show( string_SK, deliml = SK_"""" )
293 call disp%show("len(string_SK)")
294 call disp%show( len(string_SK) )
295 call disp%skip()
296
297 string_SK = SK_"AB"
298 call disp%show("string_SK")
299 call disp%show( string_SK, deliml = SK_"""" )
300 call disp%show("len(string_SK)")
301 call disp%show( len(string_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_'-')
304 call disp%show("string_SK")
305 call disp%show( string_SK, deliml = SK_"""" )
306 call disp%show("len(string_SK)")
307 call disp%show( len(string_SK) )
308 call disp%skip()
309
310 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
316 call disp%skip()
317
318 string_SK = SK_"ABC"
319 call disp%show("string_SK")
320 call disp%show( string_SK, deliml = SK_"""" )
321 call disp%show("len(string_SK)")
322 call disp%show( len(string_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_'-')
325 call disp%show("string_SK")
326 call disp%show( string_SK, deliml = SK_"""" )
327 call disp%show("len(string_SK)")
328 call disp%show( len(string_SK) )
329 call disp%skip()
330
331 string_SK = SK_"ABC"
332 call disp%show("string_SK")
333 call disp%show( string_SK, deliml = SK_"""" )
334 call disp%show("len(string_SK)")
335 call disp%show( len(string_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_'-')
338 call disp%show("string_SK")
339 call disp%show( string_SK, deliml = SK_"""" )
340 call disp%show("len(string_SK)")
341 call disp%show( len(string_SK) )
342 call disp%skip()
343
344 string_SK = SK_"AB"
345 call disp%show("string_SK")
346 call disp%show( string_SK, deliml = SK_"""" )
347 call disp%show("len(string_SK)")
348 call disp%show( len(string_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_'-')
351 call disp%show("string_SK")
352 call disp%show( string_SK, deliml = SK_"""" )
353 call disp%show("len(string_SK)")
354 call disp%show( len(string_SK) )
355 call disp%skip()
356
357 string_SK = SK_"AB"
358 call disp%show("string_SK")
359 call disp%show( string_SK, deliml = SK_"""" )
360 call disp%show("len(string_SK)")
361 call disp%show( len(string_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_'-')
364 call disp%show("string_SK")
365 call disp%show( string_SK, deliml = SK_"""" )
366 call disp%show("len(string_SK)")
367 call disp%show( len(string_SK) )
368 call disp%skip()
369
370 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
376 call disp%skip()
377
378 string_SK = SK_"In God We Trust"
379 call disp%show("string_SK")
380 call disp%show( string_SK , deliml = SK_"""" )
381 call disp%show("len(string_SK)")
382 call disp%show( len(string_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_'*')
385 call disp%show("string_SK")
386 call disp%show( string_SK, deliml = SK_"""" )
387 call disp%show("len(string_SK)")
388 call disp%show( len(string_SK) )
389 call disp%skip()
390
391 string_SK = SK_"ABCDEF"
392 call disp%show("string_SK")
393 call disp%show( string_SK, deliml = SK_"""" )
394 call disp%show("len(string_SK)")
395 call disp%show( len(string_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_'+')
398 call disp%show("string_SK")
399 call disp%show( string_SK, deliml = SK_"""" )
400 call disp%show("len(string_SK)")
401 call disp%show( len(string_SK) )
402 call disp%skip()
403
404 string_SK = SK_"ABCDEF"
405 call disp%show("string_SK")
406 call disp%show( string_SK, deliml = SK_"""" )
407 call disp%show("len(string_SK)")
408 call disp%show( len(string_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_'+')
411 call disp%show("string_SK")
412 call disp%show( string_SK, deliml = SK_"""" )
413 call disp%show("len(string_SK)")
414 call disp%show( len(string_SK) )
415 call disp%skip()
416
417contains
418
419 subroutine reset()
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)
425 end subroutine
426
427end 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 and center an array to a new larger size.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%
10! Center character scalar.
11!%%%%%%%%%%%%%%%%%%%%%%%%%
12
13string_SK
14"ABCDEF"
15len(string_SK)
16+6
17string_SK = getCentered(string_SK, size = 8_IK, fill = SK_'-')
18string_SK
19"-ABCDEF-"
20len(string_SK)
21+8
22
23!%%%%%%%%%%%%%%%%%%%%%%%%
24! Center character array.
25!%%%%%%%%%%%%%%%%%%%%%%%%
26
27array_SK
28"AA", "BB", "CC", "DD", "EE", "FF"
29lbound(array_SK), ubound(array_SK)
30+3, +8
31array_SK = getCentered(array_SK, size = 8_IK, fill = SK_'%%')
32array_SK
33"%%", "AA", "BB", "CC", "DD", "EE", "FF", "%%"
34lbound(array_SK), ubound(array_SK)
35+1, +8
36
37!%%%%%%%%%%%%%%%%%%%%%%
38! Center logical array.
39!%%%%%%%%%%%%%%%%%%%%%%
40
41array_LK
42T, T, T, T, T, T
43lbound(array_LK), ubound(array_LK)
44+3, +8
45array_LK = getCentered(array_LK, size = 8_IK, fill = .false._LK)
46array_LK
47F, T, T, T, T, T, T, F
48lbound(array_LK), ubound(array_LK)
49+1, +8
50
51!%%%%%%%%%%%%%%%%%%%%%%
52! Center integer array.
53!%%%%%%%%%%%%%%%%%%%%%%
54
55array_IK
56+1, +2, +3, +4, +5, +6
57lbound(array_IK), ubound(array_IK)
58+3, +8
59array_IK = getCentered(array_IK, size = 8_IK, fill = -9999_IK)
60array_IK
61-9999, +1, +2, +3, +4, +5, +6, -9999
62lbound(array_IK), ubound(array_IK)
63+1, +8
64
65!%%%%%%%%%%%%%%%%%%%%%%
66! Center complex array.
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)
71lbound(array_CK), ubound(array_CK)
72+3, +8
73array_CK = getCentered(array_CK, size = 8_IK, fill = (-999._CK,-999._CK))
74array_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)
77+1, +8
78
79!%%%%%%%%%%%%%%%%%%%
80! Center real array.
81!%%%%%%%%%%%%%%%%%%%
82
83array_RK
84+1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000
85lbound(array_RK), ubound(array_RK)
86+3, +8
87array_RK = getCentered(array_RK, size = 8_IK, fill = 0._RK)
88array_RK
89+0.0000000000000000, +1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000, +0.0000000000000000
90lbound(array_RK), ubound(array_RK)
91+1, +8
92
93!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
94!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
95! Shrink an array to a new smaller size. Note that fill has no effect because the array shrinks.
96!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
97!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
98
99
100!%%%%%%%%%%%%%%%%%%%%%%%%%
101! Center character scalar.
102!%%%%%%%%%%%%%%%%%%%%%%%%%
103
104string_SK
105"ABCDEF"
106len(string_SK)
107+6
108string_SK = getCentered(string_SK, size = 3_IK, fill = '~')
109string_SK
110"BCD"
111len(string_SK)
112+3
113
114!%%%%%%%%%%%%%%%%%%%%%%%%
115! Center character array.
116!%%%%%%%%%%%%%%%%%%%%%%%%
117
118array_SK
119"AA", "BB", "CC", "DD", "EE", "FF"
120lbound(array_SK), ubound(array_SK)
121+3, +8
122array_SK = getCentered(array_SK, size = 3_IK, fill = '**')
123array_SK
124"BB", "CC", "DD"
125lbound(array_SK), ubound(array_SK)
126+1, +3
127
128!%%%%%%%%%%%%%%%%%%%%%%
129! Center logical array.
130!%%%%%%%%%%%%%%%%%%%%%%
131
132array_LK
133T, T, T, T, T, T
134lbound(array_LK), ubound(array_LK)
135+3, +8
136array_LK = getCentered(array_LK, size = 3_IK, fill = .false._LK)
137array_LK
138T, T, T
139lbound(array_LK), ubound(array_LK)
140+1, +3
141
142!%%%%%%%%%%%%%%%%%%%%%%
143! Center integer array.
144!%%%%%%%%%%%%%%%%%%%%%%
145
146array_IK
147+1, +2, +3, +4, +5, +6
148lbound(array_IK), ubound(array_IK)
149+3, +8
150array_IK = getCentered(array_IK, size = 3_IK, fill = 0_IK)
151array_IK
152+2, +3, +4
153lbound(array_IK), ubound(array_IK)
154+1, +3
155
156!%%%%%%%%%%%%%%%%%%%%%%
157! Center complex array.
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)
162lbound(array_CK), ubound(array_CK)
163+3, +8
164array_CK = getCentered(array_CK, size = 3_IK, fill = (0._CK, 0._CK))
165array_CK
166(+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+4.0000000000000000, -4.0000000000000000)
167lbound(array_CK), ubound(array_CK)
168+1, +3
169
170!%%%%%%%%%%%%%%%%%%%
171! Center real array.
172!%%%%%%%%%%%%%%%%%%%
173
174array_RK
175+1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000
176lbound(array_RK), ubound(array_RK)
177+3, +8
178array_RK = getCentered(array_RK, size = 3_IK, fill = +0._RK)
179array_RK
180+2.0000000000000000, +3.0000000000000000, +4.0000000000000000
181lbound(array_RK), ubound(array_RK)
182+1, +3
183
184!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
185!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
186! Example of perfect array centering.
187!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
188!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
189
190string_SK
191"ABC"
192len(string_SK)
193+3
194string_SK = getCentered(string_SK, size = 7_IK, fill = SK_'-')
195string_SK
196"--ABC--"
197len(string_SK)
198+7
199
200string_SK
201"ABC"
202len(string_SK)
203+3
204string_SK = getCentered(string_SK, size = 1_IK, fill = SK_'-')
205string_SK
206"B"
207len(string_SK)
208+1
209
210string_SK
211"AB"
212len(string_SK)
213+2
214string_SK = getCentered(string_SK, size = 6_IK, fill = SK_'-')
215string_SK
216"--AB--"
217len(string_SK)
218+6
219
220
221!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
222!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
223! Example of imperfect array centering.
224!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
225!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
226
227string_SK
228"ABC"
229len(string_SK)
230+3
231string_SK = getCentered(string_SK, size = 6_IK, fill = SK_'-')
232string_SK
233"-ABC--"
234len(string_SK)
235+6
236
237string_SK
238"ABC"
239len(string_SK)
240+3
241string_SK = getCentered(string_SK, size = 2_IK, fill = SK_'-')
242string_SK
243"AB"
244len(string_SK)
245+2
246
247string_SK
248"AB"
249len(string_SK)
250+2
251string_SK = getCentered(string_SK, size = 3_IK, fill = SK_'-')
252string_SK
253"AB-"
254len(string_SK)
255+3
256
257string_SK
258"AB"
259len(string_SK)
260+2
261string_SK = getCentered(string_SK, size = 1_IK, fill = SK_'-')
262string_SK
263"A"
264len(string_SK)
265+1
266
267
268!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
269!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
270! Centering array contents within a given size, padded by left and right margins.
271!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
272!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
273
274string_SK
275"In God We Trust"
276len(string_SK)
277+15
278string_SK = getCentered(string_SK, size = 45_IK, fill = SK_' ', lmsize = 4_IK, rmsize = 4_IK, lmfill = SK_'*', rmfill = SK_'*')
279string_SK
280"**** In God We Trust ****"
281len(string_SK)
282+53
283
284string_SK
285"ABCDEF"
286len(string_SK)
287+6
288string_SK = getCentered(string_SK, size = 8_IK, lmsize = 3_IK, rmsize = 5_IK, fill = SK_'-', lmfill = SK_'*', rmfill = SK_'+')
289string_SK
290"***-ABCDEF-+++++"
291len(string_SK)
292+16
293
294string_SK
295"ABCDEF"
296len(string_SK)
297+6
298string_SK = getCentered(string_SK, size = 2_IK, lmsize = 3_IK, rmsize = 5_IK, fill = SK_'-', lmfill = SK_'*', rmfill = SK_'+')
299string_SK
300"***CD+++++"
301len(string_SK)
302+10
303
304
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.

  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 at Austin

Definition at line 169 of file pm_arrayCenter.F90.


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