Line data Source code
1 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3 : !!!!
4 : !!!! MIT License
5 : !!!!
6 : !!!! ParaMonte: plain powerful parallel Monte Carlo library.
7 : !!!!
8 : !!!! Copyright (C) 2012-present, The Computational Data Science Lab
9 : !!!!
10 : !!!! This file is part of the ParaMonte library.
11 : !!!!
12 : !!!! Permission is hereby granted, free of charge, to any person obtaining a
13 : !!!! copy of this software and associated documentation files (the "Software"),
14 : !!!! to deal in the Software without restriction, including without limitation
15 : !!!! the rights to use, copy, modify, merge, publish, distribute, sublicense,
16 : !!!! and/or sell copies of the Software, and to permit persons to whom the
17 : !!!! Software is furnished to do so, subject to the following conditions:
18 : !!!!
19 : !!!! The above copyright notice and this permission notice shall be
20 : !!!! included in all copies or substantial portions of the Software.
21 : !!!!
22 : !!!! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
23 : !!!! EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
24 : !!!! MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
25 : !!!! IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM,
26 : !!!! DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
27 : !!!! OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE
28 : !!!! OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
29 : !!!!
30 : !!!! ACKNOWLEDGMENT
31 : !!!!
32 : !!!! ParaMonte is an honor-ware and its currency is acknowledgment and citations.
33 : !!!! As per the ParaMonte library license agreement terms, if you use any parts of
34 : !!!! this library for any purposes, kindly acknowledge the use of ParaMonte in your
35 : !!!! work (education/research/industry/development/...) by citing the ParaMonte
36 : !!!! library as described on this page:
37 : !!!!
38 : !!!! https://github.com/cdslaborg/paramonte/blob/main/ACKNOWLEDGMENT.md
39 : !!!!
40 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
41 : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
42 :
43 : !> \brief This module contains procedures for sorting arrays.
44 : !> \author Amir Shahmoradi
45 :
46 : module Sort_mod
47 :
48 : implicit none
49 :
50 : character(*), parameter :: MODULE_NAME = "@Sort_mod"
51 :
52 : interface sortArray
53 : module procedure :: sortArray_RK, sortAscending_RK, sortAscendingWithRooter_RK
54 : end interface sortArray
55 :
56 : interface sortAscending
57 : module procedure :: sortAscending_RK, sortAscendingWithRooter_RK
58 : end interface sortAscending
59 :
60 : interface indexArray
61 : module procedure :: indexArray_IK, indexArray_RK
62 : end interface indexArray
63 :
64 : interface getMedian
65 : module procedure :: getMedian_RK
66 : end interface getMedian
67 :
68 : contains
69 :
70 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
71 :
72 : !> \brief
73 : !> Sort (recursively) the input real array of arbitrary size from the smallest value to the largest and
74 : !> return the result inside the input array.
75 : !>
76 : !> @param[inout] array : The input vector to be sorted.
77 : !>
78 : !> \warning
79 : !> On return, the contents of the input array is completely overwritten.
80 297 : pure recursive subroutine sortArray_RK(array)
81 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
82 : !DEC$ ATTRIBUTES DLLEXPORT :: sortArray_RK
83 : #endif
84 : use Constants_mod, only: IK, RK
85 : implicit none
86 : real(RK), intent(inout) :: array(:)
87 : integer(IK) :: iq
88 :
89 297 : if(size(array) > 1) then
90 147 : call partition(array, iq)
91 147 : call sortArray_RK(array(:iq-1))
92 147 : call sortArray_RK(array(iq:))
93 : endif
94 297 : end subroutine sortArray_RK
95 :
96 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
97 :
98 294 : pure subroutine partition(array, marker)
99 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
100 : !DEC$ ATTRIBUTES DLLEXPORT :: partition
101 : #endif
102 : use Constants_mod, only: IK, RK
103 : implicit none
104 : real(RK) , intent(inout) :: array(:)
105 : integer(IK), intent(out) :: marker
106 : integer(IK) :: i, j
107 147 : real(RK) :: temp
108 147 : real(RK) :: x ! pivot point
109 147 : x = array(1)
110 147 : i= 0
111 147 : j= size(array) + 1
112 225 : do
113 372 : j = j-1
114 312 : do
115 684 : if (array(j) <= x) exit
116 312 : j = j-1
117 : end do
118 372 : i = i+1
119 159 : do
120 531 : if (array(i) >= x) exit
121 159 : i = i+1
122 : end do
123 372 : if (i < j) then ! exchange array(i) and array(j)
124 225 : temp = array(i)
125 225 : array(i) = array(j)
126 225 : array(j) = temp
127 147 : elseif (i == j) then
128 39 : marker = i+1
129 39 : return
130 : else
131 108 : marker = i
132 108 : return
133 : endif
134 : end do
135 147 : end subroutine partition
136 :
137 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
138 :
139 : !> \brief
140 : !> Sort (recursively) the input real array of arbitrary size from the smallest value to the largest and
141 : !> return the result inside the input array.
142 : !>
143 : !> @param[in] np : The length of the input vector to be sorted.
144 : !> @param[inout] Point : The input vector to be sorted.
145 : !> @param[out] Err : An object of class [Err_type](@ref err_mod::err_type).
146 : !>
147 : !> \warning
148 : !> On return, the contents of the input array is completely overwritten by the output sorted array.
149 : !>
150 : !> \warning
151 : !> On return, the value of `Err%occurred` must be checked for any potential occurrences of errors during sorting.
152 633 : pure subroutine sortAscending_RK(np,Point,Err)
153 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
154 : !DEC$ ATTRIBUTES DLLEXPORT :: sortAscending_RK
155 : #endif
156 147 : use Constants_mod, only: IK, RK
157 : use Misc_mod, only: swap
158 : use Err_mod, only: Err_type
159 :
160 : implicit none
161 :
162 : integer(IK) , parameter :: nn = 15, NSTACK = 100
163 : integer(IK) , intent(in) :: np
164 : real(RK) , intent(inout) :: Point(np)
165 : type(Err_type) , intent(out) :: Err
166 :
167 : character(*) , parameter :: PROCEDURE_NAME = MODULE_NAME//"@sortAscending_RK()"
168 579 : real(RK) :: dummy
169 : integer(IK) :: k,i,j,jstack,m,r,istack(NSTACK)
170 :
171 579 : Err%occurred = .false.
172 :
173 579 : jstack=0
174 579 : m = 1
175 579 : r = np
176 140052 : do
177 140631 : if (r-m < nn) then
178 660840 : do j = m+1,r
179 590235 : dummy = Point(j)
180 1390160 : do i = j-1,m,-1
181 1349220 : if (Point(i) <= dummy) exit
182 1390160 : Point(i+1) = Point(i)
183 : end do
184 660840 : Point(i+1) = dummy
185 : end do
186 70605 : if (jstack == 0) return
187 70026 : r = istack(jstack)
188 70026 : m = istack(jstack-1)
189 70026 : jstack = jstack-2
190 : else
191 70026 : k = (m+r)/2
192 70026 : call swap(Point(k),Point(m+1))
193 70026 : if (Point(m)>Point(r)) call swap(Point(m),Point(r))
194 70026 : if (Point(m+1)>Point(r)) call swap(Point(m+1),Point(r))
195 70026 : if (Point(m)>Point(m+1)) call swap(Point(m),Point(m+1))
196 70026 : i = m+1
197 70026 : j = r
198 70026 : dummy = Point(m+1)
199 1473100 : do
200 1954930 : do
201 3498050 : i = i+1
202 3498050 : if (Point(i) >= dummy) exit
203 : end do
204 1835380 : do
205 3378500 : j = j-1
206 3378500 : if (Point(j) <= dummy) exit
207 : end do
208 1543120 : if (j < i) exit
209 1473100 : call swap(Point(i),Point(j))
210 : end do
211 70026 : Point(m+1) = Point(j)
212 70026 : Point(j) = dummy
213 70026 : jstack = jstack+2
214 70026 : if (jstack > NSTACK) then
215 : ! LCOV_EXCL_START
216 : Err%occurred = .true.
217 : Err%msg = PROCEDURE_NAME//": NSTACK is too small."
218 : return
219 : ! LCOV_EXCL_STOP
220 : end if
221 70026 : if (r-i+1 >= j-m) then
222 37221 : istack(jstack) = r
223 37221 : istack(jstack-1) = i
224 37221 : r = j-1
225 : else
226 32805 : istack(jstack) = j-1
227 32805 : istack(jstack-1) = m
228 32805 : m = i
229 : end if
230 : end if
231 : end do
232 579 : end subroutine sortAscending_RK
233 :
234 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
235 :
236 : !> \brief
237 : !> Given the real `Array(1:n)` return the array `Indx(1:n)` such that `Array(Indx(j)), j=1:n` is in ascending order.
238 : !>
239 : !> @param[in] n : The length of the input vector to be sorted.
240 : !> @param[in] Array : The input vector of length `n` to be sorted.
241 : !> @param[out] Indx : The output integer vector of indices of the sorted vector.
242 : !> @param[out] Err : An object of class [Err_type](@ref err_mod::err_type).
243 : !>
244 : !> \warning
245 : !> On return, the value of `Err%%occurred` must be checked for any potential occurrences of errors during sorting.
246 859 : pure subroutine indexArray_RK(n,Array,Indx,Err)
247 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
248 : !DEC$ ATTRIBUTES DLLEXPORT :: indexArray_RK
249 : #endif
250 579 : use Constants_mod, only: IK, RK
251 : use Misc_mod, only: swap
252 : use Err_mod, only: Err_type
253 :
254 : implicit none
255 :
256 : integer(IK) , intent(in) :: n
257 : real(RK) , intent(in) :: Array(n)
258 : integer(IK) , intent(out) :: Indx(n)
259 : type(Err_type) , intent(out) :: Err
260 :
261 : character(*) , parameter :: PROCEDURE_NAME = MODULE_NAME//"@indexArray_RK()"
262 : integer(IK) , parameter :: nn=15, NSTACK=50
263 : integer(IK) :: k,i,j,indext,jstack,l,r
264 : integer(IK) :: istack(NSTACK)
265 859 : real(RK) :: a
266 :
267 859 : Err%occurred = .false.
268 :
269 236893 : do j = 1,n
270 236893 : Indx(j) = j
271 : end do
272 859 : jstack=0
273 859 : l=1
274 859 : r=n
275 47159 : do
276 47159 : if (r-l < nn) then
277 212884 : do j=l+1,r
278 188875 : indext=Indx(j)
279 188875 : a=Array(indext)
280 689497 : do i=j-1,l,-1
281 649062 : if (Array(Indx(i)) <= a) exit
282 689497 : Indx(i+1)=Indx(i)
283 : end do
284 212884 : Indx(i+1)=indext
285 : end do
286 24009 : if (jstack == 0) return
287 23150 : r=istack(jstack)
288 23150 : l=istack(jstack-1)
289 23150 : jstack=jstack-2
290 : else
291 23150 : k=(l+r)/2
292 23150 : call swap(Indx(k),Indx(l+1))
293 23150 : call exchangeIndex(Indx(l),Indx(r))
294 23150 : call exchangeIndex(Indx(l+1),Indx(r))
295 23150 : call exchangeIndex(Indx(l),Indx(l+1))
296 23150 : i=l+1
297 23150 : j=r
298 23150 : indext=Indx(l+1)
299 23150 : a=Array(indext)
300 264568 : do
301 394554 : do
302 682272 : i=i+1
303 682272 : if (Array(Indx(i)) >= a) exit
304 : end do
305 377316 : do
306 665034 : j=j-1
307 665034 : if (Array(Indx(j)) <= a) exit
308 : end do
309 287718 : if (j < i) exit
310 264568 : call swap(Indx(i),Indx(j))
311 : end do
312 23150 : Indx(l+1)=Indx(j)
313 23150 : Indx(j)=indext
314 23150 : jstack=jstack+2
315 23150 : if (jstack > NSTACK) then
316 : ! LCOV_EXCL_START
317 : Err%occurred = .true.
318 : Err%msg = PROCEDURE_NAME//": NSTACK is too small."
319 : return
320 : ! LCOV_EXCL_STOP
321 : end if
322 23150 : if (r-i+1 >= j-l) then
323 11518 : istack(jstack)=r
324 11518 : istack(jstack-1)=i
325 11518 : r=j-1
326 : else
327 11632 : istack(jstack)=j-1
328 11632 : istack(jstack-1)=l
329 11632 : l=i
330 : end if
331 : end if
332 : end do
333 : contains
334 69450 : pure subroutine exchangeIndex(i,j)
335 : integer(IK), intent(inout) :: i,j
336 : integer(IK) :: swp
337 69450 : if (Array(j) < Array(i)) then
338 26844 : swp=i
339 26844 : i=j
340 26844 : j=swp
341 : end if
342 859 : end subroutine exchangeIndex
343 : end subroutine indexArray_RK
344 :
345 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
346 :
347 : !> \brief
348 : !> Given the integer `Array(1:n)` return the array `Indx(1:n)` such that `Array(Indx(j)), j=1:n` is in ascending order.
349 : !>
350 : !> @param[in] n : The length of the input vector to be sorted.
351 : !> @param[in] Array : The input vector of length `n` to be sorted.
352 : !> @param[out] Indx : The output integer vector of indices of the sorted vector.
353 : !> @param[out] Err : An object of class [Err_type](@ref err_mod::err_type).
354 : !>
355 : !> \warning
356 : !> On return, the value of `Err%%occurred` must be checked for any potential occurrences of errors during sorting.
357 218 : pure subroutine indexArray_IK(n,Array,Indx,Err)
358 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
359 : !DEC$ ATTRIBUTES DLLEXPORT :: indexArray_IK
360 : #endif
361 : use Constants_mod, only: IK, RK
362 : use Misc_mod, only: swap
363 : use Err_mod, only: Err_type
364 :
365 : implicit none
366 :
367 : integer(IK) , intent(in) :: n
368 : integer(IK) , intent(in) :: Array(n)
369 : integer(IK) , intent(out) :: Indx(n)
370 : type(Err_type) , intent(out) :: Err
371 :
372 : character(*) , parameter :: PROCEDURE_NAME = MODULE_NAME//"@indexArray_IK()"
373 : integer(IK) , parameter :: nn=15_IK, NSTACK=50_IK
374 : integer(IK) :: k,i,j,indext,jstack,l,r
375 : integer(IK) :: istack(NSTACK)
376 : integer(IK) :: a
377 :
378 218 : Err%occurred = .false.
379 :
380 1039 : do j = 1,n
381 1039 : Indx(j) = j
382 : end do
383 218 : jstack=0
384 218 : l=1
385 218 : r=n
386 242 : do
387 242 : if (r-l < nn) then
388 809 : do j=l+1,r
389 579 : indext=Indx(j)
390 579 : a=Array(indext)
391 991 : do i=j-1,l,-1
392 961 : if (Array(Indx(i)) <= a) exit
393 991 : Indx(i+1)=Indx(i)
394 : end do
395 809 : Indx(i+1)=indext
396 : end do
397 230 : if (jstack == 0) return
398 12 : r=istack(jstack)
399 12 : l=istack(jstack-1)
400 12 : jstack=jstack-2
401 : else
402 12 : k=(l+r)/2
403 12 : call swap(Indx(k),Indx(l+1))
404 12 : call exchangeIndex(Indx(l),Indx(r))
405 12 : call exchangeIndex(Indx(l+1),Indx(r))
406 12 : call exchangeIndex(Indx(l),Indx(l+1))
407 12 : i=l+1
408 12 : j=r
409 12 : indext=Indx(l+1)
410 12 : a=Array(indext)
411 66 : do
412 174 : do
413 252 : i=i+1
414 252 : if (Array(Indx(i)) >= a) exit
415 : end do
416 63 : do
417 141 : j=j-1
418 141 : if (Array(Indx(j)) <= a) exit
419 : end do
420 78 : if (j < i) exit
421 66 : call swap(Indx(i),Indx(j))
422 : end do
423 12 : Indx(l+1)=Indx(j)
424 12 : Indx(j)=indext
425 12 : jstack=jstack+2
426 12 : if (jstack > NSTACK) then
427 : ! LCOV_EXCL_START
428 : Err%occurred = .true.
429 : Err%msg = PROCEDURE_NAME//": NSTACK is too small."
430 : return
431 : ! LCOV_EXCL_STOP
432 : end if
433 12 : if (r-i+1 >= j-l) then
434 3 : istack(jstack)=r
435 3 : istack(jstack-1)=i
436 3 : r=j-1
437 : else
438 9 : istack(jstack)=j-1
439 9 : istack(jstack-1)=l
440 9 : l=i
441 : end if
442 : end if
443 : end do
444 : contains
445 36 : pure subroutine exchangeIndex(i,j)
446 : integer(IK), intent(inout) :: i,j
447 : integer(IK) :: swp
448 36 : if (Array(j) < Array(i)) then
449 18 : swp=i
450 18 : i=j
451 18 : j=swp
452 : end if
453 218 : end subroutine exchangeIndex
454 : end subroutine indexArray_IK
455 :
456 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
457 :
458 : !> \brief
459 : !> Sort the real `Leader(1:lenLeader)` in ascending order using Quicksort while making the corresponding rearrangement of real `Rooter(1:lenLeader)`.
460 : !>
461 : !> @param[in] lenLeader : The length of the input vector to be sorted.
462 : !> @param[inout] Leader : The vector of length `lenLeader` to be sorted.
463 : !> @param[inout] Rooter : The vector of length `lenLeader` to be sorted according to the rearrangement of the elements of `Leader`.
464 : !> @param[out] Err : An object of class [Err_type](@ref err_mod::err_type).
465 : !>
466 : !> \warning
467 : !> On return, the value of `Err%%occurred` must be checked for any potential occurrences of errors during sorting.
468 3 : pure subroutine sortAscendingWithRooter_IK(lenLeader,Leader,Rooter,Err)
469 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
470 : !DEC$ ATTRIBUTES DLLEXPORT :: sortAscendingWithRooter_IK
471 : #endif
472 : use Constants_mod, only: IK
473 : use Err_mod, only: Err_type
474 :
475 : implicit none
476 :
477 : integer(IK) , intent(in) :: lenLeader
478 : integer(IK) , intent(inout) :: Leader(lenLeader), Rooter(lenLeader)
479 : type(Err_type) , intent(out) :: Err
480 :
481 : character(*) , parameter :: PROCEDURE_NAME = MODULE_NAME//"@sortAscendingWithRooter_RK()"
482 3 : integer(IK) :: Indx(lenLeader)
483 :
484 3 : call indexArray_IK(lenLeader,Leader,Indx,Err)
485 : ! LCOV_EXCL_START
486 : if (Err%occurred) then
487 : Err%msg = PROCEDURE_NAME//": NSTACK is too small."
488 : return
489 : end if
490 : ! LCOV_EXCL_STOP
491 303 : Leader = Leader(Indx)
492 303 : Rooter = Rooter(Indx)
493 3 : end subroutine sortAscendingWithRooter_IK
494 :
495 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
496 :
497 : !> \brief
498 : !> Sort the real `Leader(1:lenLeader)` in ascending order using Quicksort while making the corresponding rearrangement of real `Rooter(1:lenLeader)`.
499 : !>
500 : !> @param[in] lenLeader : The length of the input vector to be sorted.
501 : !> @param[inout] Leader : The vector of length `lenLeader` to be sorted.
502 : !> @param[inout] Rooter : The vector of length `lenLeader` to be sorted according to the rearrangement of the elements of `Leader`.
503 : !> @param[out] Err : An object of class [Err_type](@ref err_mod::err_type).
504 : !>
505 : !> \warning
506 : !> On return, the value of `Err%%occurred` must be checked for any potential occurrences of errors during sorting.
507 15 : pure subroutine sortAscendingWithRooter_RK(lenLeader,Leader,Rooter,Err)
508 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
509 : !DEC$ ATTRIBUTES DLLEXPORT :: sortAscendingWithRooter_RK
510 : #endif
511 3 : use Constants_mod, only: IK, RK
512 : use Err_mod, only: Err_type
513 :
514 : implicit none
515 :
516 : integer(IK) , intent(in) :: lenLeader
517 : real(RK) , intent(inout) :: Leader(lenLeader), Rooter(lenLeader)
518 : type(Err_type) , intent(out) :: Err
519 :
520 : character(*) , parameter :: PROCEDURE_NAME = MODULE_NAME//"@sortAscendingWithRooter_RK()"
521 15 : integer(IK) :: Indx(lenLeader)
522 :
523 15 : call indexArray_RK(lenLeader,Leader,Indx,Err)
524 : ! LCOV_EXCL_START
525 : if (Err%occurred) then
526 : Err%msg = PROCEDURE_NAME//": NSTACK is too small."
527 : return
528 : end if
529 : ! LCOV_EXCL_STOP
530 1515 : Leader = Leader(Indx)
531 1515 : Rooter = Rooter(Indx)
532 15 : end subroutine sortAscendingWithRooter_RK
533 :
534 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
535 :
536 : !> \brief
537 : !> Return the median of the input vector.
538 : !>
539 : !> @param[in] lenArray : The length of the input vector.
540 : !> @param[in] Array : The input vector of length `lenArray` whose median is to be found.
541 : !> @param[out] median : The median of the input `Array`.
542 : !> @param[out] Err : An object of class [Err_type](@ref err_mod::err_type).
543 : !>
544 : !> \warning
545 : !> On return, the value of `Err%%occurred` must be checked for any potential occurrences of errors during sorting.
546 96 : pure subroutine getMedian_RK(lenArray,Array,median,Err)
547 : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
548 : !DEC$ ATTRIBUTES DLLEXPORT :: getMedian_RK
549 : #endif
550 15 : use Constants_mod, only: IK, RK
551 : use Err_mod, only: Err_type
552 :
553 : implicit none
554 :
555 : integer(IK) , intent(in) :: lenArray
556 : real(RK) , intent(in) :: Array(lenArray)
557 : real(RK) , intent(out) :: median
558 : type(Err_type) , intent(out) :: Err
559 :
560 : character(*) , parameter :: PROCEDURE_NAME = MODULE_NAME//"@median_RK()"
561 663 : real(RK) :: ArrayDummy(lenArray)
562 : integer(IK) :: lenArrayHalf
563 :
564 663 : ArrayDummy = Array
565 96 : call sortAscending(np=lenArray,Point=ArrayDummy,Err=Err)
566 : ! LCOV_EXCL_START
567 : if (Err%occurred) then
568 : Err%msg = PROCEDURE_NAME//Err%msg
569 : return
570 : end if
571 : ! LCOV_EXCL_STOP
572 :
573 96 : lenArrayHalf = lenArray / 2
574 96 : median = ArrayDummy(lenArrayHalf+1)
575 96 : if (2*lenArrayHalf==lenArray) median = 0.5_RK * ( ArrayDummy(lenArrayHalf) + median )
576 :
577 96 : end subroutine getMedian_RK
578 :
579 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
580 :
581 0 : end module Sort_mod
|