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

Center the contents of the input array within the output arrayCentered while filling the newly added elements (if any) with the user-specified fill. More...

Detailed Description

Center the contents of the input array within the output arrayCentered while filling the newly added elements (if any) with the user-specified fill.

Additionally, if the left and right margins lmsize, rmsize are specified, center the input array within the section arrayCentered(lmsize + 1 : size(arrayCentered) - rmsize) and optionally fill the margins with the user-specified lmfill and rmfill.

Parameters
[in,out]arrayCentered: The output contiguous array-like object of shape (:) of the same type and kind as the input array, containing the contents of array centered within it, with the requested margins and fillings.
[in,out]array: The input contiguous array of shape (:) of either
  • type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) of the same length type-parameter as that of array, 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 character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU),
whose contents will be copied to the center of arrayCentered after subtracting the lmsize and rmsize elements from the left and the right of arrayCentered.
[in]lmsize: The input scalar integer of default kind IK representing the width of the left-margin of the output arrayCentered
(optional, default = 0. It must be present if and only if the input argument rmsize is also present.)
[in]rmsize: The input scalar integer of default kind IK representing the width of the right-margin of the output arrayCentered
(optional, default = 0. It must be present if and 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 arrayCentered 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 in arrayCentered 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 output arrayCentered. 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 output arrayCentered. 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)


Possible calling interfaces

call setCentered(arrayCentered, array, fill = fill)
call setCentered(arrayCentered, array, lmsize, rmsize, fill = fill, lmfill = lmfill, rmfill = rmfill)
Center the contents of the input array within the output arrayCentered while filling the newly added ...
This module contains procedures and generic interfaces for resizing an input array and centering the ...
Warning
Note that the new elements of the output arrayCentered are not initialized to any particular value if fill is missing.
Therefore, the contents of the new elements of the output array 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 (after taking into account the left and the right margins).
In such a case, there will be no new elements to initialize their values and fill will have no effects on the output.
The sum of the input lmsize and rmsize arguments must not be larger than the size of the input arrayCentered.
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.
See also
setResized
setRebound
setRefilled
setRebilled
getCoreHalo
setCoreHalo
getCentered
setCentered
getPadded
setPadded


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

Status: Unresolved
Source: GNU Fortran Compiler gfortran version 10.3
Description: GNU Fortran Compiler gfortran cannot handle regular allocation for assumed-length allocatable character types and returns the following error.
Fortran runtime error: Integer overflow when calculating the amount of memory to allocate


Remedy (as of ParaMonte Library version 1.5): The preprocessor conditions bypass this bug.

Remedy (as of ParaMonte Library version 2.0.0): The new interface allocatable output arguments entirely, thus obviating the need to handle this gracefully.

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 Austin

Definition at line 1119 of file pm_arrayCenter.F90.


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