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

Generate a copy of the input array whose elements are reordered according to the input index array such that,
  arrayNew = array(index) or,
  arrayNew = array(index(size(index):1:-1)) if action = reverse or,
holds for the output. More...

Detailed Description

Generate a copy of the input array whose elements are reordered according to the input index array such that,
  arrayNew = array(index) or,
  arrayNew = array(index(size(index):1:-1)) if action = reverse or,
holds for the output.

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 logical of kind any supported by the processor (e.g., LK), or
  • type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64), 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),
whose elements will be reordered according to the input index.
If the input argument arrayNew is specified, then array has intent(in) and will not change.
If the input argument arrayNew is missing, then array will be reallocated and reordered on output.
The the allocatable attribute and the reallocation of array is essential for efficient fast reordering of array.
[in]index: The input contiguous array of shape (:) of type integer of default kind IK of the same size as array, containing the reordering indices.
[in]action: The input scalar object that can be
  1. the constant reverse or equivalently, an object of type reverse_type.
    Specifying this value implies that the input index must be reversed (index(size(index):1:-1)) prior to being used for remapping.
    On output, the input array will be reordered to array(index(size(index):1:-1)).
(optional, If missing, the input index will be used as is, such that array will be reordered to array(index) on output.)
Returns
arrayNew : The output array of the same type, kind, shape, and size as the input array that will contain the reordered array


Possible calling interfaces

arrayNew = getRemapped(array, index) ! `array` must be allocatable.
arrayNew = getRemapped(array, index, action) ! `array` must be allocatable.
Generate a copy of the input array whose elements are reordered according to the input index array su...
This module contains procedures and generic interfaces for remapping arrays of various types.
Warning
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.
The sizes of array, arrayNew, and index arguments must be equal.
All elements of index must be within the lower and upper bounds of array.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The primary purpose of the procedures under this generic interface is to provide a convenient bug-free method of generating a remapped copy of an array of any intrinsic type and kind, similar to the what is implicitly done by current compilers with the regular Fortran syntax.
With further compiler and language template enhancements in the future, the need for the procedures under this generic interface might be resolved in the future.
See pm_arrayRemap for the relevant benchmarks.
See also
setRemapped
setShuffled
getReversed
setReversed


Example usage

1program example
2
3 use pm_kind, only: SK ! All kinds are supported.
4 use pm_kind, only: LK ! 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
9 use pm_arrayRemap, only: getRemapped, reverse
10
11 implicit none
12
13 integer(IK) , allocatable :: index(:) ! Must be of default kind IK
14 character(:, SK), allocatable :: string1_SK , string2_SK
15 character(9, SK), allocatable :: Array1_SK(:) , Array2_SK(:) ! Can be any processor-supported kind.
16 integer(IK) , allocatable :: Array1_IK(:) , Array2_IK(:) ! Can be any processor-supported kind.
17 complex(CK) , allocatable :: Array1_CK(:) , Array2_CK(:) ! Can be any processor-supported kind.
18 real(RK) , allocatable :: Array1_RK(:) , Array2_RK(:) ! Can be any processor-supported kind.
19 logical(LK) , allocatable :: Array1_LK(:) , Array2_LK(:)
20
21 type(display_type) :: disp
22
23 disp = display_type(file = "main.out.F90")
24
25 call disp%skip()
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show("! Remap an ALLOCATABLE array IN-PLACE.")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%skip()
32
33 index = [2_IK, 1_IK, 4_IK, 3_IK, 3_IK, 2_IK]
34
35 string1_SK = "ABCDEF"
36 Array1_SK = ["AA", "BB", "CC", "DD", "EE", "FF"]
37 Array1_IK = [1_IK, 2_IK, 3_IK, 4_IK, 5_IK, 6_IK]
38 Array1_RK = [1._RK, 2._RK, 3._RK, 4._RK, 5._RK, 6._RK]
39 Array1_CK = [(1._CK, -1._CK), (2._CK, -2._CK), (3._CK, -3._CK), (4._CK, -4._CK), (5._CK, -5._CK), (6._CK, -6._CK)]
40 Array1_LK = [.false._LK, .true._LK, .false._LK, .true._LK, .false._LK, .true._LK]
41
42 string2_SK = string1_SK
43 Array2_SK = Array1_SK
44 Array2_IK = Array1_IK
45 Array2_RK = Array1_RK
46 Array2_CK = Array1_CK
47 Array2_LK = Array1_LK
48
49 call disp%skip()
50 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
51 call disp%show("! Remap character scalar.")
52 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
53 call disp%skip()
54
55 call disp%show("string2_SK")
56 call disp%show( string2_SK, deliml = SK_"""" )
57 call disp%show("index")
58 call disp%show( index )
59 call disp%show("string2_SK = getRemapped(string2_SK, index)")
60 string2_SK = getRemapped(string2_SK, index)
61 call disp%show("string2_SK")
62 call disp%show( string2_SK, deliml = SK_"""" )
63
64 call disp%skip()
65 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
66 call disp%show("! Remap character array.")
67 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
68 call disp%skip()
69
70 call disp%show("Array2_SK")
71 call disp%show( Array2_SK, deliml = SK_"""" )
72 call disp%show("index")
73 call disp%show( index )
74 call disp%show("Array2_SK = getRemapped(Array2_SK, index)")
75 Array2_SK = getRemapped(Array2_SK, index)
76 call disp%show("Array2_SK")
77 call disp%show( Array2_SK, deliml = SK_"""" )
78
79 call disp%skip()
80 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
81 call disp%show("! Remap logical array.")
82 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
83 call disp%skip()
84
85 call disp%show("Array2_LK")
86 call disp%show( Array2_LK )
87 call disp%show("index")
88 call disp%show( index )
89 call disp%show("Array2_LK = getRemapped(Array2_LK, index)")
90 Array2_LK = getRemapped(Array2_LK, index)
91 call disp%show("Array2_LK")
92 call disp%show( Array2_LK )
93
94 call disp%skip()
95 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
96 call disp%show("! Remap integer array.")
97 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
98 call disp%skip()
99
100 call disp%show("Array2_IK")
101 call disp%show( Array2_IK )
102 call disp%show("index")
103 call disp%show( index )
104 call disp%show("Array2_IK = getRemapped(Array2_IK, index)")
105 Array2_IK = getRemapped(Array2_IK, index)
106 call disp%show("Array2_IK")
107 call disp%show( Array2_IK )
108
109 call disp%skip()
110 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
111 call disp%show("! Remap complex array.")
112 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
113 call disp%skip()
114
115 call disp%show("Array2_CK")
116 call disp%show( Array2_CK )
117 call disp%show("index")
118 call disp%show( index )
119 call disp%show("Array2_CK = getRemapped(Array2_CK, index)")
120 Array2_CK = getRemapped(Array2_CK, index)
121 call disp%show("Array2_CK")
122 call disp%show( Array2_CK )
123
124 call disp%skip()
125 call disp%show("!%%%%%%%%%%%%%%%%%%")
126 call disp%show("! Remap real array.")
127 call disp%show("!%%%%%%%%%%%%%%%%%%")
128 call disp%skip()
129
130 call disp%show("Array2_RK")
131 call disp%show( Array2_RK )
132 call disp%show("index")
133 call disp%show( index )
134 call disp%show("Array2_RK = getRemapped(Array2_RK, index)")
135 Array2_RK = getRemapped(Array2_RK, index)
136 call disp%show("Array2_RK")
137 call disp%show( Array2_RK )
138
139
140 call disp%skip()
141 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
142 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
143 call disp%show("! Remap an ALLOCATABLE array IN-PLACE in backward direction.")
144 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
145 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
146 call disp%skip()
147
148
149 index = [1_IK, 2_IK, 3_IK, 4_IK, 5_IK, 6_IK]
150 string2_SK = string1_SK
151 Array2_SK = Array1_SK
152 Array2_IK = Array1_IK
153 Array2_RK = Array1_RK
154 Array2_CK = Array1_CK
155 Array2_LK = Array1_LK
156
157 call disp%skip()
158 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
159 call disp%show("! Remap character scalar.")
160 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
161 call disp%skip()
162
163 call disp%show("string2_SK")
164 call disp%show( string2_SK, deliml = SK_"""" )
165 call disp%show("index")
166 call disp%show( index )
167 call disp%show("string2_SK = getRemapped(string2_SK, index, reverse)")
168 string2_SK = getRemapped(string2_SK, index, reverse)
169 call disp%show("string2_SK")
170 call disp%show( string2_SK, deliml = SK_"""" )
171
172 call disp%skip()
173 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
174 call disp%show("! Remap character array.")
175 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
176 call disp%skip()
177
178 call disp%show("Array2_SK")
179 call disp%show( Array2_SK, deliml = SK_"""" )
180 call disp%show("index")
181 call disp%show( index )
182 call disp%show("Array2_SK = getRemapped(Array2_SK, index, reverse)")
183 Array2_SK = getRemapped(Array2_SK, index, reverse)
184 call disp%show("Array2_SK")
185 call disp%show( Array2_SK, deliml = SK_"""" )
186
187 call disp%skip()
188 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
189 call disp%show("! Remap logical array.")
190 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
191 call disp%skip()
192
193 call disp%show("Array2_LK")
194 call disp%show( Array2_LK )
195 call disp%show("index")
196 call disp%show( index )
197 call disp%show("Array2_LK = getRemapped(Array2_LK, index, reverse)")
198 Array2_LK = getRemapped(Array2_LK, index, reverse)
199 call disp%show("Array2_LK")
200 call disp%show( Array2_LK )
201
202 call disp%skip()
203 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
204 call disp%show("! Remap integer array.")
205 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
206 call disp%skip()
207
208 call disp%show("Array2_IK")
209 call disp%show( Array2_IK )
210 call disp%show("index")
211 call disp%show( index )
212 call disp%show("Array2_IK = getRemapped(Array2_IK, index, reverse)")
213 Array2_IK = getRemapped(Array2_IK, index, reverse)
214 call disp%show("Array2_IK")
215 call disp%show( Array2_IK )
216
217 call disp%skip()
218 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
219 call disp%show("! Remap complex array.")
220 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
221 call disp%skip()
222
223 call disp%show("Array2_CK")
224 call disp%show( Array2_CK )
225 call disp%show("index")
226 call disp%show( index )
227 call disp%show("Array2_CK = getRemapped(Array2_CK, index, reverse)")
228 Array2_CK = getRemapped(Array2_CK, index, reverse)
229 call disp%show("Array2_CK")
230 call disp%show( Array2_CK )
231
232 call disp%skip()
233 call disp%show("!%%%%%%%%%%%%%%%%%%")
234 call disp%show("! Remap real array.")
235 call disp%show("!%%%%%%%%%%%%%%%%%%")
236 call disp%skip()
237
238 call disp%show("Array2_RK")
239 call disp%show( Array2_RK )
240 call disp%show("index")
241 call disp%show( index )
242 call disp%show("Array2_RK = getRemapped(Array2_RK, index, reverse)")
243 Array2_RK = getRemapped(Array2_RK, index, reverse)
244 call disp%show("Array2_RK")
245 call disp%show( Array2_RK )
246
247end 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! Remap an ALLOCATABLE array IN-PLACE.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%
10! Remap character scalar.
11!%%%%%%%%%%%%%%%%%%%%%%%%
12
13string2_SK
14"ABCDEF"
15index
16+2, +1, +4, +3, +3, +2
17string2_SK = getRemapped(string2_SK, index)
18string2_SK
19"BADCCB"
20
21!%%%%%%%%%%%%%%%%%%%%%%%
22! Remap character array.
23!%%%%%%%%%%%%%%%%%%%%%%%
24
25Array2_SK
26"AA ", "BB ", "CC ", "DD ", "EE ", "FF "
27index
28+2, +1, +4, +3, +3, +2
29Array2_SK = getRemapped(Array2_SK, index)
30Array2_SK
31"BB ", "AA ", "DD ", "CC ", "CC ", "BB "
32
33!%%%%%%%%%%%%%%%%%%%%%
34! Remap logical array.
35!%%%%%%%%%%%%%%%%%%%%%
36
37Array2_LK
38F, T, F, T, F, T
39index
40+2, +1, +4, +3, +3, +2
41Array2_LK = getRemapped(Array2_LK, index)
42Array2_LK
43T, F, T, F, F, T
44
45!%%%%%%%%%%%%%%%%%%%%%
46! Remap integer array.
47!%%%%%%%%%%%%%%%%%%%%%
48
49Array2_IK
50+1, +2, +3, +4, +5, +6
51index
52+2, +1, +4, +3, +3, +2
53Array2_IK = getRemapped(Array2_IK, index)
54Array2_IK
55+2, +1, +4, +3, +3, +2
56
57!%%%%%%%%%%%%%%%%%%%%%
58! Remap complex array.
59!%%%%%%%%%%%%%%%%%%%%%
60
61Array2_CK
62(+1.0000000000000000, -1.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000)
63index
64+2, +1, +4, +3, +3, +2
65Array2_CK = getRemapped(Array2_CK, index)
66Array2_CK
67(+2.0000000000000000, -2.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+2.0000000000000000, -2.0000000000000000)
68
69!%%%%%%%%%%%%%%%%%%
70! Remap real array.
71!%%%%%%%%%%%%%%%%%%
72
73Array2_RK
74+1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000
75index
76+2, +1, +4, +3, +3, +2
77Array2_RK = getRemapped(Array2_RK, index)
78Array2_RK
79+2.0000000000000000, +1.0000000000000000, +4.0000000000000000, +3.0000000000000000, +3.0000000000000000, +2.0000000000000000
80
81!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
83! Remap an ALLOCATABLE array IN-PLACE in backward direction.
84!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
85!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
86
87
88!%%%%%%%%%%%%%%%%%%%%%%%%
89! Remap character scalar.
90!%%%%%%%%%%%%%%%%%%%%%%%%
91
92string2_SK
93"ABCDEF"
94index
95+1, +2, +3, +4, +5, +6
96string2_SK = getRemapped(string2_SK, index, reverse)
97string2_SK
98"FEDCBA"
99
100!%%%%%%%%%%%%%%%%%%%%%%%
101! Remap character array.
102!%%%%%%%%%%%%%%%%%%%%%%%
103
104Array2_SK
105"AA ", "BB ", "CC ", "DD ", "EE ", "FF "
106index
107+1, +2, +3, +4, +5, +6
108Array2_SK = getRemapped(Array2_SK, index, reverse)
109Array2_SK
110"FF ", "EE ", "DD ", "CC ", "BB ", "AA "
111
112!%%%%%%%%%%%%%%%%%%%%%
113! Remap logical array.
114!%%%%%%%%%%%%%%%%%%%%%
115
116Array2_LK
117F, T, F, T, F, T
118index
119+1, +2, +3, +4, +5, +6
120Array2_LK = getRemapped(Array2_LK, index, reverse)
121Array2_LK
122T, F, T, F, T, F
123
124!%%%%%%%%%%%%%%%%%%%%%
125! Remap integer array.
126!%%%%%%%%%%%%%%%%%%%%%
127
128Array2_IK
129+1, +2, +3, +4, +5, +6
130index
131+1, +2, +3, +4, +5, +6
132Array2_IK = getRemapped(Array2_IK, index, reverse)
133Array2_IK
134+6, +5, +4, +3, +2, +1
135
136!%%%%%%%%%%%%%%%%%%%%%
137! Remap complex array.
138!%%%%%%%%%%%%%%%%%%%%%
139
140Array2_CK
141(+1.0000000000000000, -1.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000)
142index
143+1, +2, +3, +4, +5, +6
144Array2_CK = getRemapped(Array2_CK, index, reverse)
145Array2_CK
146(+6.0000000000000000, -6.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+1.0000000000000000, -1.0000000000000000)
147
148!%%%%%%%%%%%%%%%%%%
149! Remap real array.
150!%%%%%%%%%%%%%%%%%%
151
152Array2_RK
153+1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000
154index
155+1, +2, +3, +4, +5, +6
156Array2_RK = getRemapped(Array2_RK, index, reverse)
157Array2_RK
158+6.0000000000000000, +5.0000000000000000, +4.0000000000000000, +3.0000000000000000, +2.0000000000000000, +1.0000000000000000
159
Test:
test_pm_arrayRemap
Todo:
Low Priority: This generic interface can be extended to 1D input objects.
Todo:
Critical Priority: The gfortran bugs in the implementations of this generic interface must be resolved in the future.


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 140 of file pm_arrayRemap.F90.


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