Generate and return an output array whose elements are the reversed-order elements of the input array
, such that
array = array(lenArray:1:-1)
,
where lenArray = len(array)
if array
is a scalar character, or lenArray = size(array)
is an array of rank 1
.
More...
Generate and return an output array whose elements are the reversed-order elements of the input array
, such that
array = array(lenArray:1:-1)
,
where lenArray = len(array)
if array
is a scalar character, or lenArray = size(array)
is an array of rank 1
.
- Parameters
-
[in] | array | : The input contiguous array of shape (:) of either
-
type css_pdt or,
-
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
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) of arbitrary length type parameter,
whose elements will be reversed in order.
|
- Returns
ArrayReversed
: The output contiguous
array of the same type, kind, shape, and size as the input array
that will contain the reversed array.
Possible calling interfaces ⛓
Generate and return an output array whose elements are the reversed-order elements of the input array...
This module contains procedures and generic interfaces for reversing the order of elements in arrays ...
- 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.
- See also
- setReversed
setShuffled
getRemapped
setRemapped
Example usage ⛓
13 character(:, SK),
allocatable :: string1_SK , string2_SK
14 character(
9, SK),
allocatable :: Array1_SK(:) , Array2_SK(:)
15 integer(IK) ,
allocatable :: Array1_IK(:) , Array2_IK(:)
16 complex(CK) ,
allocatable :: Array1_CK(:) , Array2_CK(:)
17 real(RK) ,
allocatable :: Array1_RK(:) , Array2_RK(:)
18 logical(LK) ,
allocatable :: Array1_LK(:) , Array2_LK(:)
20 type(display_type) :: disp
25 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show(
"! Reverse an array via a function call.")
28 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 Array1_SK
= [
"AA",
"BB",
"CC",
"DD",
"EE",
"FF"]
34 Array1_IK
= [
1_IK,
2_IK,
3_IK,
4_IK,
5_IK,
6_IK]
35 Array1_RK
= [
1._RK,
2._RK,
3._RK,
4._RK,
5._RK,
6._RK]
36 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)]
37 Array1_LK
= [
.false._LK,
.true._LK,
.false._LK,
.true._LK,
.false._LK,
.true._LK]
39 string2_SK
= string1_SK
47 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 call disp%show(
"! Reverse character scalar.")
49 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%")
53 call disp%show( string2_SK, deliml
= SK_
"""" )
54 call disp%show(
"getReversed(string2_SK)")
58 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%show(
"! Reverse character array.")
60 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
64 call disp%show( Array2_SK, deliml
= SK_
"""" )
65 call disp%show(
"getReversed(Array2_SK)")
67 call disp%show( Array2_SK, deliml
= SK_
"""" )
70 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
71 call disp%show(
"! Reverse logical array.")
72 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
77 call disp%show(
"getReversed(Array2_LK)")
81 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
82 call disp%show(
"! Reverse integer array.")
83 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
88 call disp%show(
"getReversed(Array2_IK)")
92 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
93 call disp%show(
"! Reverse complex array.")
94 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%")
99 call disp%show(
"getReversed(Array2_CK)")
103 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%")
104 call disp%show(
"! Reverse real array.")
105 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%")
110 call disp%show(
"getReversed(Array2_RK)")
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 ⛓
23"AA ",
"BB ",
"CC ",
"DD ",
"EE ",
"FF "
25"FF ",
"EE ",
"DD ",
"CC ",
"BB ",
"AA "
50(
+1.0000000000000000,
-1.0000000000000000), (
+2.0000000000000000,
-2.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000), (
+5.0000000000000000,
-5.0000000000000000), (
+6.0000000000000000,
-6.0000000000000000)
52(
+6.0000000000000000,
-6.0000000000000000), (
+5.0000000000000000,
-5.0000000000000000), (
+4.0000000000000000,
-4.0000000000000000), (
+3.0000000000000000,
-3.0000000000000000), (
+2.0000000000000000,
-2.0000000000000000), (
+1.0000000000000000,
-1.0000000000000000)
59+1.0000000000000000,
+2.0000000000000000,
+3.0000000000000000,
+4.0000000000000000,
+5.0000000000000000,
+6.0000000000000000
61+6.0000000000000000,
+5.0000000000000000,
+4.0000000000000000,
+3.0000000000000000,
+2.0000000000000000,
+1.0000000000000000
- Test:
- test_pm_arrayReverse
- Todo:
- Low Priority: This generic interface can be extended to 2D input objects.
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 116 of file pm_arrayReverse.F90.