Generate a resized copy of the input array
by padding it to the right with the requested paddings and optionally adding margins to the right of the padded array
optionally filled with the corresponding fills.
- Parameters
-
[in] | 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), 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
or,
-
a scalar assumed-length
character of default kind SK.
If the rmsize input argument is missing, then the size of the output arrayPadded is lenArray + rpsize , otherwise the length of the output arrayPadded is lenArray + rmsize + rpsize . |
[in] | rpsize | : The input scalar of type integer of default kind IK representing the number of rpfill to add to the right of the array. |
[in] | rpfill | : The input scalar of the same type and kind as the input array containing the value to fill the right padding (if any) of the output arrayPadded . If array is of type character , then
-
the equality
len(rpfill) == 1 must also hold if array is a scalar string.
-
the equality
len(rpfill) == len(array) must also hold if array is an array of strings.
|
[in] | rmsize | : The input scalar integer of default kind IK representing the size of the right-margin of the output arrayPadded
(optional, default = 0 ). |
[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(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 right-margin will remain uninitialized. It can be present only if rmsize is also present). |
- Returns
arrayPadded
: The output object of the same type, kind, and rank as the input array
with the same lower bound, whose size is lenArray + rpsize + rmsize
, whose contents are the same as the contents of array
but padded to the right, optionally with the specified right margin.
Possible calling interfaces ⛓
arrayPadded
= getPaddedr(array, rpsize, rpfill, rmsize, rmfill
= rmfill)
Generate a resized copy of the input array by padding it to the right with the requested paddings and...
This module contains procedures and generic interfaces for resizing an input array and padding them w...
- Warning
- The right margin elements will not be initialized to any particular values unless the corresponding
rmfill
argument is present.
-
The input
rpsize
and 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.
- See also
- setResized
setRebound
setRefilled
setRebilled
getCoreHalo
setCoreHalo
getCentered
setCentered
getPadded
setPadded
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(
"! Pad an array to the left with symbols.")
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(
"! Pad a character scalar.")
43 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
46 string_SK
= SK_
"In God We Trust"
48 call disp%show( string_SK, deliml
= SK_
"""" )
51 call disp%show(
"string_SK = getPaddedr(string_SK, rpsize = 5_IK, rpfill = SK_'-')")
52 string_SK
= getPaddedr(string_SK, rpsize
= 5_IK, rpfill
= SK_
'-')
54 call disp%show( string_SK, deliml
= SK_
"""" )
59 string_SK
= SK_
"In Science We Invest"
61 call disp%show( string_SK, deliml
= SK_
"""" )
64 call disp%show(
"string_SK = getPaddedr(string_SK, rpsize = 4_IK, rpfill = SK_' ', rmsize = 4_IK, rmfill = SK_'*')")
65 string_SK
= getPaddedr(string_SK, rpsize
= 4_IK, rpfill
= SK_
' ', rmsize
= 4_IK, rmfill
= SK_
'*')
67 call disp%show( string_SK, deliml
= SK_
"""" )
72 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
73 call disp%show(
"! Pad a character array.")
74 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
78 call disp%show( Array_SK, deliml
= SK_
"""" )
79 call disp%show(
"lbound(Array_SK), ubound(Array_SK)")
80 call disp%show([
lbound(Array_SK),
ubound(Array_SK)])
81 call disp%show(
"Array_SK = getPaddedr(Array_SK, rpsize = 5_IK, rpfill = SK_' ')")
82 Array_SK
= getPaddedr(Array_SK, rpsize
= 5_IK, rpfill
= SK_
' ')
84 call disp%show( Array_SK, deliml
= SK_
"""" )
85 call disp%show(
"lbound(Array_SK), ubound(Array_SK)")
86 call disp%show([
lbound(Array_SK),
ubound(Array_SK)])
89 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
90 call disp%show(
"! Pad a logical array.")
91 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
96 call disp%show(
"lbound(Array_LK), ubound(Array_LK)")
97 call disp%show([
lbound(Array_LK),
ubound(Array_LK)])
98 call disp%show(
"Array_LK = getPaddedr(Array_LK, rpsize = 5_IK, rpfill = .false._LK)")
99 Array_LK
= getPaddedr(Array_LK, rpsize
= 5_IK, rpfill
= .false._LK)
102 call disp%show(
"lbound(Array_LK), ubound(Array_LK)")
103 call disp%show([
lbound(Array_LK),
ubound(Array_LK)])
106 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
107 call disp%show(
"! Pad an integer array.")
108 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
113 call disp%show(
"lbound(Array_IK), ubound(Array_IK)")
114 call disp%show([
lbound(Array_IK),
ubound(Array_IK)])
115 call disp%show(
"Array_IK = getPaddedr(Array_IK, rpsize = 5_IK, rpfill = -0_IK)")
116 Array_IK
= getPaddedr(Array_IK, rpsize
= 5_IK, rpfill
= -0_IK)
119 call disp%show(
"lbound(Array_IK), ubound(Array_IK)")
120 call disp%show([
lbound(Array_IK),
ubound(Array_IK)])
123 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
124 call disp%show(
"! Pad a complex array.")
125 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
130 call disp%show(
"lbound(Array_CK), ubound(Array_CK)")
131 call disp%show([
lbound(Array_CK),
ubound(Array_CK)])
132 call disp%show(
"Array_CK = getPaddedr(Array_CK, rpsize = 5_IK, rpfill = (-999._CK,-999._CK))")
133 Array_CK
= getPaddedr(Array_CK, rpsize
= 5_IK, rpfill
= (
-999._CK,
-999._CK))
136 call disp%show(
"lbound(Array_CK), ubound(Array_CK)")
137 call disp%show([
lbound(Array_CK),
ubound(Array_CK)])
140 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
141 call disp%show(
"! Pad a real array.")
142 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
147 call disp%show(
"lbound(Array_RK), ubound(Array_RK)")
148 call disp%show([
lbound(Array_RK),
ubound(Array_RK)])
149 call disp%show(
"Array_RK = getPaddedr(Array_RK, rpsize = 5_IK, rpfill = -0._RK)")
150 Array_RK
= getPaddedr(Array_RK, rpsize
= 5_IK, rpfill
= -0._RK)
153 call disp%show(
"lbound(Array_RK), ubound(Array_RK)")
154 call disp%show([
lbound(Array_RK),
ubound(Array_RK)])
159 if (
allocated(Array_SK))
deallocate(Array_SK)
160 if (
allocated(Array_IK))
deallocate(Array_IK)
161 if (
allocated(Array_RK))
deallocate(Array_RK)
162 if (
allocated(Array_CK))
deallocate(Array_CK)
163 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
= getPaddedr(string_SK, rpsize
= 5_IK, rpfill
= SK_
'-')
27string_SK
= getPaddedr(string_SK, rpsize
= 4_IK, rpfill
= SK_
' ', rmsize
= 4_IK, rmfill
= SK_
'*')
29"In Science We Invest ****"
38"AA",
"BB",
"CC",
"DD",
"EE",
"FF"
39lbound(Array_SK),
ubound(Array_SK)
41Array_SK
= getPaddedr(Array_SK, rpsize
= 5_IK, rpfill
= SK_
' ')
43"AA",
"BB",
"CC",
"DD",
"EE",
"FF",
" ",
" ",
" ",
" ",
" "
44lbound(Array_SK),
ubound(Array_SK)
53lbound(Array_LK),
ubound(Array_LK)
55Array_LK
= getPaddedr(Array_LK, rpsize
= 5_IK, rpfill
= .false._LK)
57T, T, T, T, T, T, F, F, F, F, F
58lbound(Array_LK),
ubound(Array_LK)
67lbound(Array_IK),
ubound(Array_IK)
69Array_IK
= getPaddedr(Array_IK, rpsize
= 5_IK, rpfill
= -0_IK)
71+1,
+2,
+3,
+4,
+5,
+6,
+0,
+0,
+0,
+0,
+0
72lbound(Array_IK),
ubound(Array_IK)
80(
+1.0000000000000000,
-1.0000000000000000), (
+2.0000000000000000,
-2.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000), (
+5.0000000000000000,
-5.0000000000000000), (
+6.0000000000000000,
-6.0000000000000000)
81lbound(Array_CK),
ubound(Array_CK)
83Array_CK
= getPaddedr(Array_CK, rpsize
= 5_IK, rpfill
= (
-999._CK,
-999._CK))
85(
+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), (
-999.00000000000000,
-999.00000000000000), (
-999.00000000000000,
-999.00000000000000), (
-999.00000000000000,
-999.00000000000000), (
-999.00000000000000,
-999.00000000000000)
86lbound(Array_CK),
ubound(Array_CK)
94+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
95lbound(Array_RK),
ubound(Array_RK)
97Array_RK
= getPaddedr(Array_RK, rpsize
= 5_IK, rpfill
= -0._RK)
99+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000,
-0.0000000000000000,
-0.0000000000000000,
-0.0000000000000000,
-0.0000000000000000,
-0.0000000000000000
100lbound(Array_RK),
ubound(Array_RK)
- Test:
- test_pm_arrayPad
- 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 padded 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 3952 of file pm_arrayPad.F90.