Skip to content
Snippets Groups Projects
Commit aff65eb4 authored by WAUTELET Philippe's avatar WAUTELET Philippe
Browse files

Philippe 17/01/2020: move quicksort to tools + add kpos optional dummy argument

parent 4a84802a
No related branches found
No related tags found
No related merge requests found
!MNH_LIC Copyright 1995-2019 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC Copyright 1995-2020 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
......@@ -242,6 +242,7 @@ END MODULE MODI_RAIN_ICE_RED
! P. Wautelet 10/04/2019: replace ABORT and STOP calls by Print_msg
! P. Wautelet 28/05/2019: move COUNTJV function to tools.f90
! P. Wautelet 29/05/2019: remove PACK/UNPACK intrinsics (to get more performance and better OpenACC support)
! P. Wautelet 17/01/2020: move Quicksort to tools.f90
!
!* 0. DECLARATIONS
! ------------
......@@ -1729,30 +1730,5 @@ CONTAINS
!
!
END SUBROUTINE CORRECT_NEGATIVITIES
!
recursive subroutine quicksort(a, first, last)
implicit none
integer a(*), x, t
integer first, last
integer i, j
x = a( (first+last) / 2 )
i = first
j = last
do
do while (a(i) < x)
i=i+1
end do
do while (x < a(j))
j=j-1
end do
if (i >= j) exit
t = a(i); a(i) = a(j); a(j) = t
i=i+1
j=j-1
end do
if (first < i-1) call quicksort(a, first, i-1)
if (j+1 < last) call quicksort(a, j+1, last)
end subroutine quicksort
END SUBROUTINE RAIN_ICE_RED
!MNH_LIC Copyright 2019-2019 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC Copyright 2019-2020 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
......@@ -19,12 +19,14 @@ module mode_tools
!
! Modifications:
! P. Wautelet 28/05/2019: move COUNTJV function to tools.f90
! P. Wautelet 17/01/2020: move Quicksort to tools.f90
implicit none
private
public :: Countjv
public :: Quicksort
public :: Upcase
interface Countjv
......@@ -79,6 +81,44 @@ function Countjv3d(ltab,i1,i2,i3) result(ic)
end function Countjv3d
recursive subroutine Quicksort( ka, kbeg, kend, kpos )
integer, dimension(:), intent(inout) :: ka
integer, intent(in) :: kbeg, kend
integer, dimension(:), optional, intent(inout) :: kpos
integer :: ji, jj
integer :: itmp, itmp2, ival
ival = ka( ( kbeg + kend ) / 2 )
ji = kbeg
jj = kend
do
do while ( ka(ji) < ival )
ji = ji + 1
end do
do while ( ival < ka(jj) )
jj = jj - 1
end do
if ( ji >= jj ) exit
itmp = ka(ji)
ka(ji) = ka(jj)
ka(jj) = itmp
if ( present( kpos ) ) then
itmp2 = kpos(ji)
kpos(ji) = kpos(jj)
kpos(jj) = itmp2
end if
ji=ji+1
jj=jj-1
end do
if ( kbeg < ji - 1 ) call Quicksort( ka, kbeg, ji - 1, kpos )
if ( jj + 1 < kend ) call Quicksort( ka, jj + 1, kend, kpos )
end subroutine Quicksort
function Upcase(hstring)
character(len=*), intent(in) :: hstring
character(len=len(hstring)) :: upcase
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment