Reorder the elements of the input array
according to the input index
array, such that
array = array(index)
or,
array = array(index(size(index):1:-1))
if action =
reverse or,
arrayNew = array(index)
if arrayNew
is specified as input argument or,
arrayNew = array(index(size(index):1:-1))
if arrayNew
and action =
reverse are specified as input arguments,
- Parameters
-
[in,out] | array | : The input contiguous or input/output 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),
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
-
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.) |
[out] | arrayNew | : The output contiguous array of the same type, kind, shape, and size as the input array that will contain the reordered array.
(optional, if missing, the ordered array will stored and output in the input allocatable array in-place) |
Possible calling interfaces ⛓
!
Reorder the elements of the input array according to the input index array, such that   array = ar...
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
.
- Note
- Upon return, if
array
is an allocatable
with intent(inout)
, then it is guaranteed to have the same lower bound as before. This happens when the output argument arrayNew
is missing.
- See also
- getRemapped
setShuffled
getReversed
setReversed
Example usage ⛓
13 integer(IK) ,
allocatable ::
index(:)
14 character(:, SK),
allocatable :: string1_SK , string2_SK
15 character(
9, SK),
allocatable :: Array1_SK(:) , Array2_SK(:)
16 integer(IK) ,
allocatable :: Array1_IK(:) , Array2_IK(:)
17 complex(CK) ,
allocatable :: Array1_CK(:) , Array2_CK(:)
18 real(RK) ,
allocatable :: Array1_RK(:) , Array2_RK(:)
19 logical(LK) ,
allocatable :: Array1_LK(:) , Array2_LK(:)
21 type(display_type) :: disp
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 index
= [
2_IK,
1_IK,
4_IK,
3_IK,
3_IK,
2_IK]
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]
42 string2_SK
= string1_SK
50 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
51 call disp%show(
"! Remap character scalar.")
52 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
56 call disp%show( string2_SK, deliml
= SK_
"""" )
59 call disp%show(
"call setRemapped(string2_SK, index)")
62 call disp%show( string2_SK, deliml
= SK_
"""" )
65 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
66 call disp%show(
"! Remap character array.")
67 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
71 call disp%show( Array2_SK, deliml
= SK_
"""" )
74 call disp%show(
"call setRemapped(Array2_SK, index)")
77 call disp%show( Array2_SK, deliml
= SK_
"""" )
80 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
81 call disp%show(
"! Remap logical array.")
82 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
89 call disp%show(
"call setRemapped(Array2_LK, index)")
95 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
96 call disp%show(
"! Remap integer array.")
97 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
104 call disp%show(
"call setRemapped(Array2_IK, index)")
110 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
111 call disp%show(
"! Remap complex array.")
112 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
119 call disp%show(
"call setRemapped(Array2_CK, index)")
125 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
126 call disp%show(
"! Remap real array.")
127 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
134 call disp%show(
"call setRemapped(Array2_RK, index)")
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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
158 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
159 call disp%show(
"! Remap character scalar.")
160 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
164 call disp%show( string2_SK, deliml
= SK_
"""" )
167 call disp%show(
"call setRemapped(string2_SK, index, reverse)")
170 call disp%show( string2_SK, deliml
= SK_
"""" )
173 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
174 call disp%show(
"! Remap character array.")
175 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
179 call disp%show( Array2_SK, deliml
= SK_
"""" )
182 call disp%show(
"call setRemapped(Array2_SK, index, reverse)")
185 call disp%show( Array2_SK, deliml
= SK_
"""" )
188 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
189 call disp%show(
"! Remap logical array.")
190 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
197 call disp%show(
"call setRemapped(Array2_LK, index, reverse)")
203 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
204 call disp%show(
"! Remap integer array.")
205 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
212 call disp%show(
"call setRemapped(Array2_IK, index, reverse)")
218 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
219 call disp%show(
"! Remap complex array.")
220 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
227 call disp%show(
"call setRemapped(Array2_CK, index, reverse)")
233 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
234 call disp%show(
"! Remap real array.")
235 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
242 call disp%show(
"call setRemapped(Array2_RK, index, reverse)")
249 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
250 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
251 call disp%show(
"! Remap input array and return it in a new array.")
252 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
253 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
258 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
259 call disp%show(
"! Remap character scalar.")
260 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%")
264 call disp%show( string1_SK, deliml
= SK_
"""" )
267 call disp%show(
"call setRemapped(string1_SK, index, reverse, arrayNew = string2_SK)")
268 call setRemapped(string1_SK, index, reverse, arrayNew
= string2_SK)
270 call disp%show( string1_SK, deliml
= SK_
"""" )
273 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
274 call disp%show(
"! Remap character array.")
275 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
279 call disp%show( Array1_SK, deliml
= SK_
"""" )
282 call disp%show(
"call setRemapped(Array1_SK, index, reverse, arrayNew = Array2_SK)")
283 call setRemapped(Array1_SK, index, reverse, arrayNew
= Array2_SK)
285 call disp%show( Array2_SK, deliml
= SK_
"""" )
288 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
289 call disp%show(
"! Remap logical array.")
290 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
297 call disp%show(
"call setRemapped(Array1_LK, index, reverse, arrayNew = Array2_LK)")
298 call setRemapped(Array1_LK, index, reverse, arrayNew
= Array2_LK)
303 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
304 call disp%show(
"! Remap integer array.")
305 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
312 call disp%show(
"call setRemapped(Array1_IK, index, reverse, arrayNew = Array2_IK)")
313 call setRemapped(Array1_IK, index, reverse, arrayNew
= Array2_IK)
318 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
319 call disp%show(
"! Remap complex array.")
320 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%")
327 call disp%show(
"call setRemapped(Array1_CK, index, reverse, arrayNew = Array2_CK)")
328 call setRemapped(Array1_CK, index, reverse, arrayNew
= Array2_CK)
333 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
334 call disp%show(
"! Remap real array.")
335 call disp%show(
"!%%%%%%%%%%%%%%%%%%")
342 call disp%show(
"call setRemapped(Array1_RK, index, reverse, arrayNew = Array2_RK)")
343 call setRemapped(Array1_RK, index, reverse, arrayNew
= Array2_RK)
349 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
350 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
351 call disp%show(
"! Remapping preserves the lower bound of the input allocatable array.")
352 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
353 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
356 index
= [
2_IK,
1_IK,
4_IK,
3_IK,
3_IK,
2_IK]
- 4
357 deallocate(Array1_IK);
358 allocate(Array1_IK(
-5 :
-5 + size(index)
- 1))
359 Array1_IK(:)
= [
1_IK,
2_IK,
3_IK,
4_IK,
5_IK,
6_IK]
361 call disp%show(
"[lbound(Array1_IK,1), ubound(Array1_IK,1)]")
362 call disp%show( [
lbound(Array1_IK,
1),
ubound(Array1_IK,
1)] )
367 call disp%show(
"call setRemapped(Array1_IK, index)")
371 call disp%show(
"[lbound(Array1_IK,1), ubound(Array1_IK,1)]")
372 call disp%show( [
lbound(Array1_IK,
1),
ubound(Array1_IK,
1)] )
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 ⛓
26"AA ",
"BB ",
"CC ",
"DD ",
"EE ",
"FF "
31"BB ",
"AA ",
"DD ",
"CC ",
"CC ",
"BB "
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)
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)
74+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
79+2.0000000000000000,
+1.0000000000000000,
+4.0000000000000000,
+3.0000000000000000,
+3.0000000000000000,
+2.0000000000000000
105"AA ",
"BB ",
"CC ",
"DD ",
"EE ",
"FF "
107+1,
+2,
+3,
+4,
+5,
+6
110"FF ",
"EE ",
"DD ",
"CC ",
"BB ",
"AA "
119+1,
+2,
+3,
+4,
+5,
+6
129+1,
+2,
+3,
+4,
+5,
+6
131+1,
+2,
+3,
+4,
+5,
+6
134+6,
+5,
+4,
+3,
+2,
+1
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)
143+1,
+2,
+3,
+4,
+5,
+6
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)
153+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
155+1,
+2,
+3,
+4,
+5,
+6
158+6.0000000000000000,
+5.0000000000000000,
+4.0000000000000000,
+3.0000000000000000,
+2.0000000000000000,
+1.0000000000000000
174+1,
+2,
+3,
+4,
+5,
+6
175call setRemapped(string1_SK, index, reverse, arrayNew
= string2_SK)
184"AA ",
"BB ",
"CC ",
"DD ",
"EE ",
"FF "
186+1,
+2,
+3,
+4,
+5,
+6
187call setRemapped(Array1_SK, index, reverse, arrayNew
= Array2_SK)
189"FF ",
"EE ",
"DD ",
"CC ",
"BB ",
"AA "
198+1,
+2,
+3,
+4,
+5,
+6
199call setRemapped(Array1_LK, index, reverse, arrayNew
= Array2_LK)
208+1,
+2,
+3,
+4,
+5,
+6
210+1,
+2,
+3,
+4,
+5,
+6
211call setRemapped(Array1_IK, index, reverse, arrayNew
= Array2_IK)
213+6,
+5,
+4,
+3,
+2,
+1
220(
+1.0000000000000000,
-1.0000000000000000), (
+2.0000000000000000,
-2.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000), (
+5.0000000000000000,
-5.0000000000000000), (
+6.0000000000000000,
-6.0000000000000000)
222+1,
+2,
+3,
+4,
+5,
+6
223call setRemapped(Array1_CK, index, reverse, arrayNew
= Array2_CK)
225(
+6.0000000000000000,
-6.0000000000000000), (
+5.0000000000000000,
-5.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+2.0000000000000000,
-2.0000000000000000), (
+1.0000000000000000,
-1.0000000000000000)
232+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
234+1,
+2,
+3,
+4,
+5,
+6
235call setRemapped(Array1_RK, index, reverse, arrayNew
= Array2_RK)
237+6.0000000000000000,
+5.0000000000000000,
+4.0000000000000000,
+3.0000000000000000,
+2.0000000000000000,
+1.0000000000000000
245[
lbound(Array1_IK,
1),
ubound(Array1_IK,
1)]
248+1,
+2,
+3,
+4,
+5,
+6
250-2,
-3,
+0,
-1,
-1,
-2
253+4,
+3,
+6,
+5,
+5,
+4
254[
lbound(Array1_IK,
1),
ubound(Array1_IK,
1)]
- Test:
- test_pm_arrayRemap
- Bug:
Status: Unresolved
Source: GNU Fortran Compiler gfortran
version 10.3
Description: There is still a bug in gfortran as of version 10.3, where the compiler fails to properly deallocate the array of origin in a call to intrinsic move_alloc
:
call move_alloc(from = arrayNew, to = array)
Remedy (as of ParaMonte Library version 2.0.0): Currently, a preprocessor fence is added specifically for gfortran to explicitly deallocate arrayNew
after the call to move_alloc()
.
This is going to yield an extra tiny performance penalty. Upon the necessary corrections to gfortran, this extra preprocessing fence must be removed.
- Bug:
Status: Unresolved
Source: GNU Fortran Compiler gfortran
version 10.3-12
Description: GNU Fortran Compiler gfortran
currently does not accept deferred-length allocatable characters as actual argument to procedures with assumed-length allocatable dummy arguments.
Remedy (as of ParaMonte Library version 2.0.0): Avoid passing deferred-length allocatable
scalar character
arguments in the interfaces of the library until this GNU Fortran Compiler gfortran
bug is resolved.
- Bug:
Status: Unresolved, possibly Resolved in Intel Classic Fortran Compiler ifort
version 2023.
Source: Intel Classic Fortran Compiler ifort
version 2021.4
Description: The Intel Classic Fortran Compiler ifort
cannot compile interfaces with contiguous arguments whose lower bounds depend on the lower bounds of allocatable
input arguments.
The following is an illustration of the bug,
PURE module lower_bound
interface lower
PURE module subroutine test_lower(array, indices)
integer, intent(inout), allocatable :: array(:)
integer, intent(in), contiguous :: indices(lbound(array,1):)
end subroutine test_lower
end interface
end module lower_bound
submodule (lower_bound) smod
contains
PURE module procedure test_lower
end procedure
end submodule smod
Compile with
ifort lower_bound.F90 lower_bound@smod.F90 -c
to reproduce the compiler bug.
Remedy (as of ParaMonte Library version 2.0.0): For now, avoid such interfaces in the library.
- 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.
-
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 1071 of file pm_arrayRemap.F90.