@@ -11,6 +11,7 @@ submodule (stdlib_stats) stdlib_stats_median
11
11
! that are already partly sorted. While it is slightly slower for random arrays,
12
12
! ord_sort seems a better overall choice.
13
13
use stdlib_sorting, only: sort => ord_sort
14
+ use stdlib_selection, only: select
14
15
implicit none
15
16
16
17
contains
@@ -24,6 +25,7 @@ contains
24
25
real(${o1}$) :: res
25
26
26
27
integer(kind = int64) :: c, n
28
+ ${t1}$ :: val, val1
27
29
${t1}$, allocatable :: x_tmp(:)
28
30
29
31
if (.not.optval(mask, .true.) .or. size(x) == 0) then
@@ -43,16 +45,18 @@ contains
43
45
44
46
x_tmp = reshape(x, [n])
45
47
46
- call sort (x_tmp)
48
+ call select (x_tmp, c, val )
47
49
48
50
if (mod(n, 2_int64) == 0) then
51
+ call select(x_tmp, c+1, val1)
49
52
#:if t1[0] == 'r'
50
- res = sum(x_tmp(c:c+1) ) / 2._${o1}$
53
+ res = (val + val1 ) / 2._${o1}$
51
54
#:else
52
- res = sum( real(x_tmp(c:c+1), kind=${o1}$) ) / 2._${o1}$
55
+ res = (real(val, kind=${o1}$) + &
56
+ real(val1, kind=${o1}$)) / 2._${o1}$
53
57
#:endif
54
58
else
55
- res = x_tmp(c)
59
+ res = val
56
60
end if
57
61
58
62
end function ${name}$
@@ -74,6 +78,7 @@ contains
74
78
integer :: j${fj}$
75
79
#:endfor
76
80
#:endif
81
+ ${t1}$ :: val, val1
77
82
${t1}$, allocatable :: x_tmp(:)
78
83
79
84
if (.not.optval(mask, .true.) .or. size(x) == 0) then
@@ -107,17 +112,18 @@ contains
107
112
end if
108
113
#:endif
109
114
110
- call sort (x_tmp)
115
+ call select (x_tmp, c, val )
111
116
112
117
if (mod(n, 2) == 0) then
118
+ call select(x_tmp, c+1, val1)
113
119
res${reduce_subvector('j', rank, fi)}$ = &
114
120
#:if t1[0] == 'r'
115
- sum(x_tmp(c:c+1) ) / 2._${o1}$
121
+ (val + val1 ) / 2._${o1}$
116
122
#:else
117
- sum (real(x_tmp(c:c+1) , kind=${o1}$) ) / 2._${o1}$
123
+ (real(val , kind=${o1}$) + real(val1, kind=${o1}$) ) / 2._${o1}$
118
124
#:endif
119
125
else
120
- res${reduce_subvector('j', rank, fi)}$ = x_tmp(c)
126
+ res${reduce_subvector('j', rank, fi)}$ = val
121
127
end if
122
128
#:for fj in range(1, rank)
123
129
end do
@@ -141,6 +147,7 @@ contains
141
147
real(${o1}$) :: res
142
148
143
149
integer(kind = int64) :: c, n
150
+ ${t1}$ :: val, val1
144
151
${t1}$, allocatable :: x_tmp(:)
145
152
146
153
if (any(shape(x) .ne. shape(mask))) then
@@ -156,21 +163,26 @@ contains
156
163
157
164
x_tmp = pack(x, mask)
158
165
159
- call sort(x_tmp)
160
-
161
166
n = size(x_tmp, kind=int64)
162
- c = floor( (n + 1) / 2._${o1}$, kind=int64)
163
167
164
168
if (n == 0) then
165
169
res = ieee_value(1._${o1}$, ieee_quiet_nan)
166
- else if (mod(n, 2_int64) == 0) then
170
+ return
171
+ end if
172
+
173
+ c = floor( (n + 1) / 2._${o1}$, kind=int64)
174
+
175
+ call select(x_tmp, c, val)
176
+
177
+ if (mod(n, 2_int64) == 0) then
178
+ call select(x_tmp, c+1, val1)
167
179
#:if t1[0] == 'r'
168
- res = sum(x_tmp(c:c+1) ) / 2._${o1}$
180
+ res = (val + val1 ) / 2._${o1}$
169
181
#:else
170
- res = sum (real(x_tmp(c:c+1) , kind=${o1}$)) / 2._${o1}$
182
+ res = (real(val, kind=${o1}$) + real(val1 , kind=${o1}$)) / 2._${o1}$
171
183
#:endif
172
184
else if (mod(n, 2_int64) == 1) then
173
- res = x_tmp(c)
185
+ res = val
174
186
end if
175
187
176
188
end function ${name}$
@@ -192,6 +204,7 @@ contains
192
204
integer :: j${fj}$
193
205
#:endfor
194
206
#:endif
207
+ ${t1}$ :: val, val1
195
208
${t1}$, allocatable :: x_tmp(:)
196
209
197
210
if (any(shape(x) .ne. shape(mask))) then
@@ -220,23 +233,28 @@ contains
220
233
end if
221
234
#:endif
222
235
223
- call sort(x_tmp)
224
-
225
236
n = size(x_tmp, kind=int64)
226
- c = floor( (n + 1) / 2._${o1}$, kind=int64 )
227
237
228
238
if (n == 0) then
229
239
res${reduce_subvector('j', rank, fi)}$ = &
230
240
ieee_value(1._${o1}$, ieee_quiet_nan)
231
- else if (mod(n, 2_int64) == 0) then
241
+ return
242
+ end if
243
+
244
+ c = floor( (n + 1) / 2._${o1}$, kind=int64 )
245
+
246
+ call select(x_tmp, c, val)
247
+
248
+ if (mod(n, 2_int64) == 0) then
249
+ call select(x_tmp, c+1, val1)
232
250
res${reduce_subvector('j', rank, fi)}$ = &
233
251
#:if t1[0] == 'r'
234
- sum(x_tmp(c:c+1) ) / 2._${o1}$
252
+ (val + val1 ) / 2._${o1}$
235
253
#:else
236
- sum (real(x_tmp(c:c+1) , kind=${o1}$)) / 2._${o1}$
254
+ (real(val, kind=${o1}$) + real(val1 , kind=${o1}$)) / 2._${o1}$
237
255
#:endif
238
256
else if (mod(n, 2_int64) == 1) then
239
- res${reduce_subvector('j', rank, fi)}$ = x_tmp(c)
257
+ res${reduce_subvector('j', rank, fi)}$ = val
240
258
end if
241
259
242
260
deallocate(x_tmp)
0 commit comments