The ParaMonte Documentation Website
Current view: top level - kernel - Sort_mod.f90 (source / functions) Hit Total Coverage
Test: ParaMonte 1.5.1 :: Serial Kernel - Code Coverage Report Lines: 213 214 99.5 %
Date: 2021-01-08 13:03:42 Functions: 10 10 100.0 %
Legend: Lines: hit not hit

          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          99 :     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          99 :         if(size(array) > 1) then
      90          49 :             call partition(array, iq)
      91          49 :             call sortArray_RK(array(:iq-1))
      92          49 :             call sortArray_RK(array(iq:))
      93             :         endif
      94          99 :     end subroutine sortArray_RK
      95             : 
      96             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      97             : 
      98          98 :     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          49 :         real(RK)                   :: temp
     108          49 :         real(RK)                   :: x      ! pivot point
     109          49 :         x = array(1)
     110          49 :         i= 0
     111          49 :         j= size(array) + 1
     112          75 :         do
     113         124 :             j = j-1
     114         104 :             do
     115         228 :                 if (array(j) <= x) exit
     116         104 :                 j = j-1
     117             :             end do
     118         124 :             i = i+1
     119          53 :             do
     120         177 :                 if (array(i) >= x) exit
     121          53 :                 i = i+1
     122             :             end do
     123         124 :             if (i < j) then ! exchange array(i) and array(j)
     124          75 :                 temp = array(i)
     125          75 :                 array(i) = array(j)
     126          75 :                 array(j) = temp
     127          49 :             elseif (i == j) then
     128          13 :                 marker = i+1
     129          13 :                 return
     130             :             else
     131          36 :                 marker = i
     132          36 :                 return
     133             :             endif
     134             :         end do
     135          49 :   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          48 :     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          49 :         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          48 :         real(RK)                        :: dummy
     169             :         integer(IK)                     :: k,i,j,jstack,m,r,istack(NSTACK)
     170             : 
     171          48 :         Err%occurred = .false.
     172             : 
     173          48 :         jstack=0
     174          48 :         m = 1
     175          48 :         r = np
     176          60 :         do
     177         108 :             if (r-m < nn) then
     178         442 :                 do j = m+1,r
     179         364 :                     dummy = Point(j)
     180        1067 :                     do i = j-1,m,-1
     181         965 :                         if (Point(i) <=  dummy) exit
     182         805 :                         Point(i+1) = Point(i)
     183             :                     end do
     184         442 :                     Point(i+1) = dummy
     185             :                 end do
     186          78 :                 if (jstack == 0) return
     187          30 :                 r = istack(jstack)
     188          30 :                 m = istack(jstack-1)
     189          30 :                 jstack = jstack-2
     190             :             else
     191          30 :                 k = (m+r)/2
     192          30 :                 call swap(Point(k),Point(m+1))
     193          30 :                 if (Point(m)>Point(r)) call swap(Point(m),Point(r))
     194          30 :                 if (Point(m+1)>Point(r)) call swap(Point(m+1),Point(r))
     195          30 :                 if (Point(m)>Point(m+1)) call swap(Point(m),Point(m+1))
     196          30 :                 i = m+1
     197          30 :                 j = r
     198          30 :                 dummy = Point(m+1)
     199         139 :                 do
     200         214 :                     do
     201         383 :                         i = i+1
     202         383 :                         if (Point(i) >= dummy) exit
     203             :                     end do
     204         330 :                     do
     205         499 :                         j = j-1
     206         499 :                         if (Point(j) <= dummy) exit
     207             :                     end do
     208         169 :                     if (j < i) exit
     209         139 :                     call swap(Point(i),Point(j))
     210             :                 end do
     211          30 :                 Point(m+1) = Point(j)
     212          30 :                 Point(j) = dummy
     213          30 :                 jstack = jstack+2
     214          30 :                 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          30 :                 if (r-i+1 >= j-m) then
     222          19 :                     istack(jstack) = r
     223          19 :                     istack(jstack-1) = i
     224          19 :                     r = j-1
     225             :                 else
     226          11 :                     istack(jstack) = j-1
     227          11 :                     istack(jstack-1) = m
     228          11 :                     m = i
     229             :                 end if
     230             :             end if
     231             :         end do
     232          48 :     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         709 :     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          48 :         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         709 :         real(RK)                        :: a
     266             : 
     267         709 :         Err%occurred = .false.
     268             : 
     269      177630 :         do j = 1,n
     270      177630 :             Indx(j) = j
     271             :         end do
     272         709 :         jstack=0
     273         709 :         l=1
     274         709 :         r=n
     275       34813 :         do
     276       34813 :             if (r-l < nn) then
     277      159869 :                 do j=l+1,r
     278      142108 :                     indext=Indx(j)
     279      142108 :                     a=Array(indext)
     280      524331 :                     do i=j-1,l,-1
     281      493919 :                         if (Array(Indx(i)) <= a) exit
     282      412635 :                         Indx(i+1)=Indx(i)
     283             :                     end do
     284      159869 :                     Indx(i+1)=indext
     285             :                 end do
     286       17761 :                 if (jstack == 0) return
     287       17052 :                 r=istack(jstack)
     288       17052 :                 l=istack(jstack-1)
     289       17052 :                 jstack=jstack-2
     290             :             else
     291       17052 :                 k=(l+r)/2
     292       17052 :                 call swap(Indx(k),Indx(l+1))
     293       17052 :                 call exchangeIndex(Indx(l),Indx(r))
     294       17052 :                 call exchangeIndex(Indx(l+1),Indx(r))
     295       17052 :                 call exchangeIndex(Indx(l),Indx(l+1))
     296       17052 :                 i=l+1
     297       17052 :                 j=r
     298       17052 :                 indext=Indx(l+1)
     299       17052 :                 a=Array(indext)
     300      188343 :                 do
     301      284195 :                     do
     302      489590 :                         i=i+1
     303      489590 :                         if (Array(Indx(i)) >= a) exit
     304             :                     end do
     305      266086 :                     do
     306      471481 :                         j=j-1
     307      471481 :                         if (Array(Indx(j)) <= a) exit
     308             :                     end do
     309      205395 :                     if (j < i) exit
     310      188343 :                     call swap(Indx(i),Indx(j))
     311             :                 end do
     312       17052 :                 Indx(l+1)=Indx(j)
     313       17052 :                 Indx(j)=indext
     314       17052 :                 jstack=jstack+2
     315       17052 :                 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       17052 :                 if (r-i+1 >= j-l) then
     323        8501 :                     istack(jstack)=r
     324        8501 :                     istack(jstack-1)=i
     325        8501 :                     r=j-1
     326             :                 else
     327        8551 :                     istack(jstack)=j-1
     328        8551 :                     istack(jstack-1)=l
     329        8551 :                     l=i
     330             :                 end if
     331             :             end if
     332             :         end do
     333             :     contains
     334       51156 :         pure subroutine exchangeIndex(i,j)
     335             :             integer(IK), intent(inout) :: i,j
     336             :             integer(IK)                :: swp
     337       51156 :             if (Array(j) < Array(i)) then
     338       19887 :                 swp=i
     339       19887 :                 i=j
     340       19887 :                 j=swp
     341             :             end if
     342         709 :         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           3 :     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           3 :         Err%occurred = .false.
     379             : 
     380          69 :         do j = 1,n
     381          69 :             Indx(j) = j
     382             :         end do
     383           3 :         jstack=0
     384           3 :         l=1
     385           3 :         r=n
     386          11 :         do
     387          11 :             if (r-l < nn) then
     388          62 :                 do j=l+1,r
     389          55 :                     indext=Indx(j)
     390          55 :                     a=Array(indext)
     391         168 :                     do i=j-1,l,-1
     392         158 :                         if (Array(Indx(i)) <= a) exit
     393         123 :                         Indx(i+1)=Indx(i)
     394             :                     end do
     395          62 :                     Indx(i+1)=indext
     396             :                 end do
     397           7 :                 if (jstack == 0) return
     398           4 :                 r=istack(jstack)
     399           4 :                 l=istack(jstack-1)
     400           4 :                 jstack=jstack-2
     401             :             else
     402           4 :                 k=(l+r)/2
     403           4 :                 call swap(Indx(k),Indx(l+1))
     404           4 :                 call exchangeIndex(Indx(l),Indx(r))
     405           4 :                 call exchangeIndex(Indx(l+1),Indx(r))
     406           4 :                 call exchangeIndex(Indx(l),Indx(l+1))
     407           4 :                 i=l+1
     408           4 :                 j=r
     409           4 :                 indext=Indx(l+1)
     410           4 :                 a=Array(indext)
     411          22 :                 do
     412          58 :                     do
     413          84 :                         i=i+1
     414          84 :                         if (Array(Indx(i)) >= a) exit
     415             :                     end do
     416          21 :                     do
     417          47 :                         j=j-1
     418          47 :                         if (Array(Indx(j)) <= a) exit
     419             :                     end do
     420          26 :                     if (j < i) exit
     421          22 :                     call swap(Indx(i),Indx(j))
     422             :                 end do
     423           4 :                 Indx(l+1)=Indx(j)
     424           4 :                 Indx(j)=indext
     425           4 :                 jstack=jstack+2
     426           4 :                 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           4 :                 if (r-i+1 >= j-l) then
     434           1 :                     istack(jstack)=r
     435           1 :                     istack(jstack-1)=i
     436           1 :                     r=j-1
     437             :                 else
     438           3 :                     istack(jstack)=j-1
     439           3 :                     istack(jstack-1)=l
     440           3 :                     l=i
     441             :                 end if
     442             :             end if
     443             :         end do
     444             :     contains
     445          12 :         pure subroutine exchangeIndex(i,j)
     446             :             integer(IK), intent(inout) :: i,j
     447             :             integer(IK)                :: swp
     448          12 :             if (Array(j) < Array(i)) then
     449           6 :                 swp=i
     450           6 :                 i=j
     451           6 :                 j=swp
     452             :             end if
     453           3 :         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           1 :     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           1 :         integer(IK)                     :: Indx(lenLeader)
     483             : 
     484           1 :         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         101 :         Leader = Leader(Indx)
     492         101 :         Rooter = Rooter(Indx)
     493           1 :     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           5 :     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           1 :         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           5 :         integer(IK)                     :: Indx(lenLeader)
     522             : 
     523           5 :         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         505 :         Leader = Leader(Indx)
     531         505 :         Rooter = Rooter(Indx)
     532           5 :     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          43 :     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           5 :         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         265 :         real(RK)                        :: ArrayDummy(lenArray)
     562             :         integer(IK)                     :: lenArrayHalf
     563             : 
     564         265 :         ArrayDummy = Array
     565          43 :         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          43 :         lenArrayHalf = lenArray / 2
     574          43 :         median = ArrayDummy(lenArrayHalf+1)
     575          43 :         if (2*lenArrayHalf==lenArray) median = 0.5_RK * ( ArrayDummy(lenArrayHalf) + median )
     576             : 
     577          43 :     end subroutine getMedian_RK
     578             : 
     579             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     580             : 
     581           0 : end module Sort_mod

ParaMonte: Plain Powerful Parallel Monte Carlo Library 
The Computational Data Science Lab
© Copyright 2012 - 2021