valnth Interface

public interface valnth

NAME

orderval(3f) - [M_orderpack:FRACTILE] Return VALUE of Nth ordered
               element of array (Quick-Sort-like)

SYNOPSIS

 Function OrderVal (INVALS, NORD)

  ${TYPE} (Kind=${KIND}), Intent (In) :: INVALS(:)
  ${TYPE} (Kind=${KIND})              :: orderval
  Integer, Intent (In)                :: NORD

Where ${TYPE}(kind=${KIND}) may be

   o Real(kind=real32)
   o Real(kind=real64)
   o Integer(kind=int32)

DESCRIPTION

ORDERVAL(3f) returns the NORDth (ascending order) value of INVALS, i.e. the fractile of order NORD/SIZE(INVALS).

Internally, this subroutine simply calls ORDERLOC(3f).

This routine uses a pivoting strategy such as the one of finding the median based on the Quick-Sort algorithm, but we skew the pivot choice to try to bring it to NORD as fast as possible. It uses two temporary arrays, where it stores the indices of the values smaller than the pivot (ILOWT), and the indices of values larger than the pivot that we might still need later on (IHIGT). It iterates until it can bring the number of values in ILOWT to exactly NORD, and then finds the maximum of this set.

OPTIONS

 INVALS    array to search
 NORD     Nth lowest value to find

RETURNS

 orderval   Nth lowest value

EXAMPLES

Sample program:

program demo_orderval
!  Return value of Nth lowest value of array
use M_orderpack, only : orderval
implicit none
character(len=*),parameter :: list= '(*(g0:,", "))'
character(len=*),parameter :: sp='(*(g0,1x))'
real,parameter ::  INVALS(*)=[1.1,20.20,3.3,10.10,5.5,4.4,2.2]
integer :: i
integer :: imiddle
   write(*,list) 'ORIGINAL:',INVALS
   ! can return the same values as intrinsics minval(3f) and maxval(3f)
   print sp, 'minval',orderval(INVALS,1),          minval(INVALS)
   print sp, 'maxval',orderval(INVALS,size(INVALS)), maxval(INVALS)
   ! but more generally it can return the Nth lowest value.
   print sp,'nord=',4, ' fractile=',orderval(INVALS,4)
   ! so a value at the middle would be
   imiddle=(size(INVALS)+1)/2
   print sp,'median',orderval(INVALS,imiddle)
   ! sorting the hard way
   do i=1,size(INVALS)
      write(*,list)i,orderval(INVALS,i)
   enddo
end program demo_orderval

Results:

ORIGINAL:, 1.1000, 20.200, 3.300, 10.100, 5.500, 4.400, 2.200
minval 1.100 1.100
maxval 20.200 20.200
nord= 4  fractile= 4.400
median 4.400
1, 1.100
2, 2.200
3, 3.300
4, 4.400
5, 5.500
6, 10.100
7, 20.200

AUTHOR

Michel Olagnon - Aug. 2000

MAINTAINER

John Urban, 2022.04.16

LICENSE

CC0-1.0

Contents


Module Procedures

private function real64_valnth(INVALS, NORD) result(valnth)

!!! CASE DEFAULT !!! write (unit=,fmt=) “Assertion failed” !!! STOP

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in), Dimension (:) :: INVALS
integer, intent(in) :: NORD

Return Value real(kind=real64)

private function real32_valnth(INVALS, NORD) result(valnth)

!!! CASE DEFAULT !!! write (unit=,fmt=) “Assertion failed” !!! STOP

Arguments

Type IntentOptional Attributes Name
real(kind=real32), intent(in), Dimension (:) :: INVALS
integer, intent(in) :: NORD

Return Value real(kind=real32)

private function int32_valnth(INVALS, NORD) result(valnth)

!!! CASE DEFAULT !!! write (unit=,fmt=) “Assertion failed” !!! STOP

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in), Dimension (:) :: INVALS
integer, intent(in) :: NORD

Return Value integer(kind=int32)