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

Generate and return the Fractional rank of the input scalar string or contiguous array of rank 1 in ascending order or in the order specified by the input procedure isSorted() using the Quicksort algorithm such that array(rank) will be in ascending order (or in the requested order as specified by isSorted(). More...

Detailed Description

Generate and return the Fractional rank of the input scalar string or contiguous array of rank 1 in ascending order or in the order specified by the input procedure isSorted() using the Quicksort algorithm such that array(rank) will be in ascending order (or in the requested order as specified by isSorted().

This kind of ranking of values is widely known as fractional (1 2.5 2.5 4) ranking.
In Fractional ranking, items that compare equal receive the same ranking number, which is the mean of what they would have under ordinal rankings; equivalently, the ranking number of 1 plus the number of items ranked above it plus half the number of items equal to it.
This strategy has the property that the sum of the ranking numbers is the same as under ordinal ranking.
For this reason, it is used in computing Borda counts and ranking statistics (e.g., Spearman Correlation).
Thus if A ranks ahead of B and C (which compare equal) which are both ranked ahead of D, then A gets ranking number 1 (first), B and C each get ranking number 2.5 (average of joint second/third) and D gets ranking number 4 (fourth).
That is, if A < B == C < D, then the sequence ABCD has the Fractional ranking 1223.
Example:
Suppose the data set is 1.0, 1.0, 2.0, 3.0, 3.0, 4.0, 5.0, 5.0, 5.0.
The ordinal ranks are 1, 2, 3, 4, 5, 6, 7, 8, 9.
For v = 1.0, the Fractional rank is the average of the ordinal ranks: (1 + 2) / 2 = 1.5.
In a similar manner, for v = 5.0, the Fractional rank is (7 + 8 + 9) / 3 = 8.0.
Thus the Fractional ranks are: 1.5, 1.5, 3.0, 4.5, 4.5, 6.0, 8.0, 8.0, 8.0

Parameters
[in]array: The input contiguous array of rank 1 of either
  1. type css_pdt (parameterized container of string of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU)) or,
  2. type css_type (container of string of default kind SK) or,
  3. type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) of arbitrary length type parameter or,
  4. type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64) or,
  5. type logical of kind any supported by the processor (e.g., LK) or,
  6. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128) or,
  7. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
or,
  1. a scalar of type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) of arbitrary length type parameter,
whose elements rankings will be computed and returned.
isSorted: The external user-specified function that takes two input scalar arguments of the same type and kind as the input array.
It returns a scalar logical of default kind LK that is .true. if the first input scalar argument is sorted with respect to the second input argument according to the user-defined sorting condition within isSorted(), otherwise, it is .false..
If array is a Fortran string (i.e., a scalar character), then both input arguments to isSorted() are single character(1,SKG) where SKG is the kind of array.
The following illustrates the generic interface of isSorted() when the rank of the input array is 1,
function isSorted(a,b) result (sorted)
use pm_kind, only: SK, IK, LK, CK, RK
TYPE(KIND) , intent(in) :: a, b
logical(LK) :: sorted
end function
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
where TYPE(KIND) represents the type and kind of the input argument array, which can be one of the following,
use pm_kind, only: SK, IK, LK, CK, RK
character(*, SK) , intent(in) :: a, b
integer(IK) , intent(in) :: a, b
logical(LK) , intent(in) :: a, b
complex(CK) , intent(in) :: a, b
real(RK) , intent(in) :: a, b
type(css_type) , intent(in) :: a, b
type(css_pdt(SK)) , intent(in) :: a, b
This module contains the derived types for generating allocatable containers of scalar,...
This is the css_pdt parameterized type for generating instances of container of scalar of string obje...
This is the css_type type for generating instances of container of scalar of string objects.
where the kinds SK, IK, LK, CK, RK, can refer to any kind type parameter that is supported by the processor.
The following illustrates the generic interface of isSorted() when the input array is a scalar string,
function isSorted(a,b) result (sorted)
character(1,SKG), intent(in) :: a, b
logical(LK) :: sorted
end function
where SKG represents the kind of the input string argument array.
This user-defined equivalence check is extremely useful where a user-defined sorting criterion other than simple ascending order is needed, for example, when the case-sensitivity of an input string or array of strings is irrelevant or when sorting of the absolute values matters excluding the signs of the numbers, or when descending order is desired.
In such cases, user can define a custom sorting condition within the user-defined external function isSorted to achieve the goal.
(optional, the default sorting condition is ascending order, that is a < b.)
Returns
rank(1:size(array) : The output contiguous array of rank 1 of type real of default kind RK containing the ranks of the corresponding elements of array.
The size of rank must match that of array (or its length type parameter if array is a scalar string).
Read rank(i) as the Fractional rank of the ith element of array.


Possible calling interfaces

rank(1:size(array)) = getRankFractional(array)
rank(1:size(array)) = getRankFractional(array, isSorted)
Generate and return the Fractional rank of the input scalar string or contiguous array of rank 1 in a...
This module contains procedures and generic interfaces for obtaining various rankings of elements of ...
Warning
Note that the definition of isSorted(), if present, must be such that isSorted() .and. .not. isSorted() is equivalent to an equality check for two elements of the input array. This equality check is used to identify ties within the Standard ranking of the input array.
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 procedures under this generic interface are always impure when the input argument isSorted is present.
See also
setSelected
getRankDense
setRankDense
getRankOrdinal
setRankOrdinal
getRankFractional
setRankFractional
getRankStandard
setRankStandard
getRankModified
setRankModified
setSorted
setSorted


Example usage

1program example
2
3 use pm_io, only: display_type
4 use pm_kind, only: SK, IK, LK, CK, RK ! all other processor kinds are also supported.
5 use pm_distUnif, only: setUnifRand
7 use pm_str, only: getTrimmedTZ
8 use pm_val2str, only: getStr
9
10 implicit none
11
12 integer(IK) , parameter :: NP = 5_IK
13
14 ! Vector of indices of the sorted array.
15
16 integer(IK) :: i
17 real(RK) , allocatable :: rank(:)
18 character(:, SK) , allocatable :: string_SK
19 character(9, SK) , allocatable :: vector_SK(:)
20 real(RK) :: vector_RK(NP)
21 integer(IK) :: vector_IK(NP)
22 logical(LK) :: vector_LK(NP)
23
24 type(display_type) :: disp
25 disp = display_type(file = "main.out.F90")
26
27 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
28 ! Define the unsorted arrays.
29 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
30
31 !call setUnifRand(vector_SK, SK_"aaaaaaaaa", SK_"zzzzzzzzz")
32 call setUnifRand(vector_RK, -0.5_RK, +0.5_RK)
33 call setUnifRand(vector_IK, 1_IK, 2_IK**NP)
34 call setUnifRand(vector_LK)
35
36 call disp%skip()
37 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
38 call disp%show("!Sort string of characters in ascending order.")
39 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
40 call disp%skip()
41
42 string_SK = SK_"ParaMonte"
43 call allocateRank(len(string_SK))
44
45 call disp%skip()
46 call disp%show("string_SK")
47 call disp%show( string_SK, deliml = SK_"""" )
48 call disp%show("rank = getRankFractional(string_SK)")
49 rank = getRankFractional(string_SK)
50 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
51 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
52 call disp%skip()
53
54 call disp%skip()
55 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
56 call disp%show("!Sort array of strings of the same length in ascending order.")
57 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
58 call disp%skip()
59
60 vector_SK = [ "ParaMonte" &
61 , "V.2 " &
62 , "is " &
63 , "a " &
64 , "Parallel " &
65 , "Monte " &
66 , "Carlo " &
67 , "and " &
68 , "a " &
69 , "Machine " &
70 , "Learning " &
71 , "Library. " &
72 ]
73
74 call allocateRank(size(vector_SK))
75 call disp%skip()
76 call disp%show("vector_SK")
77 call disp%show( vector_SK, deliml = SK_"""" )
78 call disp%show("rank = getRankFractional(vector_SK)")
79 rank = getRankFractional(vector_SK)
80 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
81 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
82 call disp%skip()
83
84 call disp%skip()
85 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
86 call disp%show("!Sort vector of integers of arbitrary kinds in ascending order.")
87 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
88 call disp%skip()
89
90 call allocateRank(size(vector_IK))
91 call disp%skip()
92 call disp%show("vector_IK")
93 call disp%show( vector_IK )
94 call disp%show("rank = getRankFractional(vector_IK)")
95 rank = getRankFractional(vector_IK)
96 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
97 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
98 call disp%skip()
99
100 call allocateRank(size(vector_IK))
101 call disp%skip()
102 call disp%show("vector_IK = [1_IK, 2_IK, 3_IK, 2_IK, 1_IK]")
103 vector_IK = [1_IK, 2_IK, 3_IK, 2_IK, 1_IK]
104 call disp%show("rank = getRankFractional(vector_IK)")
105 rank = getRankFractional(vector_IK)
106 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
107 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
108 call disp%skip()
109
110 call disp%skip()
111 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
112 call disp%show("!Sort vector of logicals of arbitrary kinds in ascending order.")
113 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
114 call disp%skip()
115
116 call allocateRank(size(vector_LK))
117 call disp%skip()
118 call disp%show("vector_LK")
119 call disp%show( vector_LK )
120 call disp%show("rank = getRankFractional(vector_LK)")
121 rank = getRankFractional(vector_LK)
122 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
123 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
124 call disp%skip()
125
126 call disp%skip()
127 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
128 call disp%show("!Sort arrays of reals of arbitrary kinds in ascending order.")
129 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
130 call disp%skip()
131
132 call allocateRank(size(vector_RK))
133 call disp%skip()
134 call disp%show("vector_RK")
135 call disp%show( vector_RK )
136 call disp%show("rank = getRankFractional(vector_RK)")
137 rank = getRankFractional(vector_RK)
138 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
139 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
140 call disp%skip()
141
142 block
143 real, allocatable :: vector_RK(:)
144 vector_RK = [1.0, 1.0, 2.0, 3.0, 3.0, 4.0, 5.0, 5.0, 5.0]
145 call allocateRank(size(vector_RK))
146 call disp%skip()
147 call disp%show("vector_RK")
148 call disp%show( vector_RK )
149 call disp%show("rank = getRankFractional(vector_RK)")
150 rank = getRankFractional(vector_RK)
151 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
152 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
153 call disp%skip()
154 end block
155
156 call disp%skip()
157 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
158 call disp%show("!Sort according to an input user-defined comparison function.")
159 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
160 call disp%skip()
161
162 vector_IK = int([(i, i = 1_IK, NP)], kind = IK)
163 call allocateRank(size(vector_IK))
164 call disp%skip()
165 call disp%show("!Sort in DESCENDING (decreasing) order via an input custom-designed `isSorted()` function.")
166 call disp%skip()
167 call disp%show("vector_IK")
168 call disp%show( vector_IK )
169 call disp%show("rank = getRankFractional(vector_IK, isSorted_IK)")
170 rank = getRankFractional(vector_IK, isSorted_IK)
171 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
172 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
173 call disp%skip()
174
175 call random_number(vector_RK); vector_RK = vector_RK - 0.5_RK
176 call allocateRank(size(vector_RK))
177 call disp%skip()
178 call disp%show("!Sort in ascending order solely based on the magnitude of numbers using a custom comparison function.")
179 call disp%skip()
180 call disp%show("vector_RK")
181 call disp%show( vector_RK )
182 call disp%show("rank = getRankFractional(vector_RK, isSorted_RK)")
183 rank = getRankFractional(vector_RK, isSorted_RK)
184 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
185 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
186 call disp%skip()
187
188 string_SK = "ParaMonte"
189 call allocateRank(len(string_SK))
190 call disp%skip()
191 call disp%show("!Sort string in ascending order without case-sensitivity via a custom-designed input comparison function.")
192 call disp%skip()
193 call disp%show("string_SK")
194 call disp%show( string_SK, deliml = SK_"""" )
195 call disp%show("rank = getRankFractional(string_SK, isSorted_SK)")
196 rank = getRankFractional(string_SK, isSorted_SK)
197 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
198 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
199 call disp%skip()
200
201#if PDT_ENABLED
202 call disp%skip()
203 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
204 call disp%show("!Sort array of strings of varying length in ascending order.")
205 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
206 call disp%skip()
207
208 block
209 use pm_container, only: css_pdt
210 type(css_pdt) , allocatable :: cssvec(:)
211 cssvec =[ css_pdt("ParaMonte") &
212 , css_pdt("V.2") &
213 , css_pdt("is") &
214 , css_pdt("a") &
215 , css_pdt("Parallel") &
216 , css_pdt("Monte") &
217 , css_pdt("Carlo") &
218 , css_pdt("and") &
219 , css_pdt("a") &
220 , css_pdt("Machine") &
221 , css_pdt("Learning") &
222 , css_pdt("Library.") &
223 ]
224 call allocateRank(size(cssvec))
225 call disp%skip()
226 call disp%show("cssvec")
227 call disp%show( cssvec, deliml = SK_"""" )
228 call disp%show("rank = getRankFractional(cssvec)")
229 rank = getRankFractional(cssvec)
230 call disp%show("getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))")
231 call disp%show( getTrimmedTZ(getStr(rank, format = SK_"(*(g0,:,', '))")) )
232 call disp%skip()
233 end block
234#endif
235
236contains
237
238 function isSorted_IK(a,b) result(isSorted)
239 use pm_kind, only: LK, IK
240 integer(IK) , intent(in) :: a, b
241 logical(LK) :: isSorted
242 isSorted = a > b
243 end function
244
245 function isSorted_RK(a,b) result(isSorted)
246 use pm_kind, only: LK, IK
247 integer(IK) , intent(in) :: a, b
248 logical(LK) :: isSorted
249 isSorted = abs(a) < abs(b)
250 end function
251
252 function isSorted_SK(a,b) result(isSorted)
253 use pm_strASCII, only: getStrLower
254 use pm_kind, only: LK, SK
255 character(1, SK), intent(in) :: a, b
256 logical(LK) :: isSorted
257 isSorted = getStrLower(a) < getStrLower(b)
258 end function
259
260 subroutine allocateRank(size)
261 integer, intent(in) :: size
262 if (allocated(rank)) deallocate(rank)
263 allocate(rank(size))
264 end subroutine
265
266end program example
Return a uniform random scalar or contiguous array of arbitrary rank of randomly uniformly distribute...
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
Generate and return the input string where the uppercase English alphabets are all converted to lower...
Generate and return a vector of single-characters each element of which corresponds to one character ...
Definition: pm_str.F90:1449
Generate and return the conversion of the input value to an output Fortran string,...
Definition: pm_val2str.F90:167
This module contains classes and procedures for computing various statistical quantities related to t...
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 contains the uncommon and hardly representable ASCII characters as well as procedures for...
Definition: pm_strASCII.F90:61
This module contains classes and procedures for various string manipulations and inquiries.
Definition: pm_str.F90:49
This module contains the generic procedures for converting values of different types and kinds to For...
Definition: pm_val2str.F90:58
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!Sort string of characters in ascending order.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7string_SK
8"ParaMonte"
9rank = getRankFractional(string_SK)
10getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
112., 3.5, 8., 3.5, 1., 7., 6., 9., 5.
12
13
14!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
15!Sort array of strings of the same length in ascending order.
16!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
17
18
19vector_SK
20"ParaMonte", "V.2 ", "is ", "a ", "Parallel ", "Monte ", "Carlo ", "and ", "a ", "Machine ", "Learning ", "Library. "
21rank = getRankFractional(vector_SK)
22getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
236., 8., 12., 9.5, 7., 5., 1., 11., 9.5, 4., 2., 3.
24
25
26!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
27!Sort vector of integers of arbitrary kinds in ascending order.
28!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
29
30
31vector_IK
32+6, +32, +28, +29, +17
33rank = getRankFractional(vector_IK)
34getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
351., 5., 3., 4., 2.
36
37
38vector_IK = [1_IK, 2_IK, 3_IK, 2_IK, 1_IK]
39rank = getRankFractional(vector_IK)
40getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
411.5, 3.5, 5., 3.5, 1.5
42
43
44!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
45!Sort vector of logicals of arbitrary kinds in ascending order.
46!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
47
48
49vector_LK
50F, F, T, T, T
51rank = getRankFractional(vector_LK)
52getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
531.5, 1.5, 4., 4., 4.
54
55
56!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
57!Sort arrays of reals of arbitrary kinds in ascending order.
58!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59
60
61vector_RK
62+0.18057849512577595, -0.18096351796249610, +0.16591992473318706, -0.38846759183403767, +0.26438055969664975E-1
63rank = getRankFractional(vector_RK)
64getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
655., 2., 4., 1., 3.
66
67
68vector_RK
69+1.00000000, +1.00000000, +2.00000000, +3.00000000, +3.00000000, +4.00000000, +5.00000000, +5.00000000, +5.00000000
70rank = getRankFractional(vector_RK)
71getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
721.5, 1.5, 3., 4.5, 4.5, 6., 8., 8., 8.
73
74
75!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
76!Sort according to an input user-defined comparison function.
77!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
78
79
80!Sort in DESCENDING (decreasing) order via an input custom-designed `isSorted()` function.
81
82vector_IK
83+1, +2, +3, +4, +5
84rank = getRankFractional(vector_IK, isSorted_IK)
85getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
865., 4., 3., 2., 1.
87
88
89!Sort in ascending order solely based on the magnitude of numbers using a custom comparison function.
90
91vector_RK
92-0.17083296690461525, +0.28656530278302694, +0.10973993140036642, +0.48530256689881779, +0.20698253293536628
93rank = getRankFractional(vector_RK, isSorted_RK)
94getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
952., 1., 4., 5., 3.
96
97
98!Sort string in ascending order without case-sensitivity via a custom-designed input comparison function.
99
100string_SK
101"ParaMonte"
102rank = getRankFractional(string_SK, isSorted_SK)
103getTrimmedTZ(getStr(rank, format = SK_'(*(g0,:,', '))'))
1047., 1.5, 8., 1.5, 4., 6., 5., 9., 3.
105
106
Test:
test_pm_arrayRank
Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.5
Description: See pm_arraySplit for the description of a relevant bug in PDT name aliasing when compiled with Intel ifort 2021.5 that also applies to this module.
Remedy (as of ParaMonte Library version 2.0.0): See pm_arraySplit for the remedy.
Todo:
Low Priority: The current bypass for the PDT name aliasing bug can be reverted back to PDT name aliasing once the ifort bug is resolved.
Todo:
Low Priority: A test should be implemented for arrays of size that can be represented only by an IKD integer.


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, April 21, 2017, 1:54 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 2706 of file pm_arrayRank.F90.


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