C Library Functions  - orderval_special (3)

NAME

orderval_special(3f) - [M_orderpack:FRACTILE] Return VALUE of Nth ordered element of array (InsertSort-like)

CONTENTS

Synopsis
Description
Options
Returns
Examples
See Also
Author
Maintainer
License

SYNOPSIS

Function Orderval_Special (INVALS, INORD)

     ${TYPE} (Kind=${KIND}), Intent (In) :: INVALS(:)
     Integer, Intent (In)                :: INORD
     ${TYPE} (Kind=${KIND})              :: orderval_special

Where ${TYPE}(kind=${KIND}) may be
o Real(kind=real32)
o Real(kind=real64)
o Integer(kind=int32)

DESCRIPTION

ORDERVAL_SPECIAL(3f) returns the INORDth lowest value of INVALS(), i.e. the fractile of order INORD/SIZE(INVALS).

Internally, This subroutine uses an insertion sort, limiting insertion to the first INORD values and even less when one can know that the value that is considered will not be the INORDth.

An insertion sort is very fast when INORD is very small (2-5). Additionally, internally it requires only a work array of size INORD (and type of INVALS),

But worst case behavior can happen fairly probably (e.g., initially inverse sorted). Therefore, in many cases, the refined Quick-Sort method is faster.

so ORDERVAL_SPECIAL(3f) should be used when INORD is small and INVALS is likely to be a random array, otherwise consider using ORDERLOC(3f) or ORDERVAL(3f).

OPTIONS

INVALS input array of values
INORD specify Nth value of sorted INVALS array to return, from 1 to size(INVALS).

RETURNS

ORDERVAL_SPECIAL
  returned value

EXAMPLES

Sample program:

   program demo_orderval_special
   ! return Nth ordered value of an array
   use M_orderpack, only : orderval_special, medianval
   implicit none
   character(len=*),parameter :: list= ’(*(g0:,", "))’,sp=’(*(g0,1x))’
   integer,allocatable :: iarr(:)
   integer :: i
      iarr=[80,70,30,40,-50,60,20,10]
      print sp, ’ORIGINAL:’,iarr
      ! can return the same values as intrinsics minval(3f) and maxval(3f)
      print sp, ’minval’,orderval_special(iarr,1),          minval(iarr)
      print sp, ’maxval’,orderval_special(iarr,size(iarr)), maxval(iarr)
      ! but more generally it can return the Nth lowest value.
      print sp, ’median’,orderval_special(iarr,(size(iarr+1))/2), &
      & medianval(iarr)
      ! so only Nth ordered value can be found
      print sp,’inord=’,3, ’ fractile=’,orderval_special(iarr,3)
      ! sorting the hard way
      print sp, ’ORIGINAL:’,iarr
      do i=1,size(iarr)
         write(*,list)i,orderval_special(iarr,i)
      enddo
      print *
   end program demo_orderval_special

Results:

   ORIGINAL: 80 70 30 40 -50 60 20 10
   minval -50 -50
   maxval 80 80
   median 30 30
   inord= 3  fractile= 20
   ORIGINAL: 80 70 30 40 -50 60 20 10
   1, -50
   2, 10
   3, 20
   4, 30
   5, 40
   6, 60
   7, 70
   8, 80

SEE ALSO

ORDERLOC(3f), ORDERVAL(3f)

AUTHOR

Michel Olagnon - Aug. 2000

MAINTAINER

John Urban, 2022.04.16

LICENSE

CC0-1.0


Nemo Release 3.1 orderval_special (3) July 22, 2023
Generated by manServer 1.08 from d8e81ce3-1c99-4487-9d4f-0740948fe355 using man macros.