-
Notifications
You must be signed in to change notification settings - Fork 6
/
mo_percentile.f90
1080 lines (873 loc) · 31.9 KB
/
mo_percentile.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
MODULE mo_percentile
! This module provides routines for median and percentiles.
! Written Matthias Cuntz, Mar 2011
! Modified, Juliane Mai, Jul 2012 - different interpolation schemes in percentiles
! Matthias Cuntz, Juliane Mai, Jul 2012 - uses previous of ksmallest to half execution time
! Matthias Cuntz, May 2014 - removed numerical recipes
! License
! -------
! This file is part of the JAMS Fortran package, distributed under the MIT License.
!
! Copyright (c) 2011-2014 Matthias Cuntz, Juliane Mai, Stephan Thober - mc (at) macu (dot) de
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
USE mo_kind, ONLY: i4, sp, dp
Implicit NONE
PUBLIC :: median ! Median
PUBLIC :: n_element ! The n-th smallest value in an array
PUBLIC :: percentile ! The value below which a certain percent of the input fall
PUBLIC :: qmedian ! Quick median calculation, rearranges input
! ------------------------------------------------------------------
! NAME
! median
! PURPOSE
! Returns the median of the values in an array.
! If size is even, then the mean of the size/2 and size/2+1 element is the median.
!
! If an optinal mask is given, values only on those locations that correspond
! to true values in the mask are used.
! CALLING SEQUENCE
! out = median(vec,mask=mask)
! INTENT(IN)
! real(sp/dp) :: vec(:) 1D-array with input numbers
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp) :: out median of values in input array
! INTENT(IN), OPTIONAL
! logical :: mask(:) 1D-array of logical values with size(vec).
! If present, only those locations in vec corresponding to the true values in mask are used.
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! vec = (/ 1.,2.,3.,4.,5.,6.,7.,8.,9.,10. /)
! ! Returns 5.5
! out = median(vec)
! -> see also example in test directory
! LITERATURE
! None
! HISTORY
! Written, Matthias Cuntz, Mar 2011
! Modified, Matthias Cuntz, Juliane Mai, Jul 2012 - uses previous of ksmallest to half execution time
INTERFACE median
MODULE PROCEDURE median_sp, median_dp
END INTERFACE median
! ------------------------------------------------------------------
! NAME
! n_element
! PURPOSE
! Returns the n-th smallest value in an array.
!
! If an optinal mask is given, values only on those locations that correspond
! to true values in the mask are used.
! CALLING SEQUENCE
! out = n_element(vec, n, mask=mask, before=before, previous=previous, after=after, next=next)
! INTENT(IN)
! real(sp/dp) :: vec(:) 1D-array with input numbers
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp) :: out n-th smallest value in input array
! INTENT(IN), OPTIONAL
! integer(i4) :: n index of sorted array
! logical :: mask(:) 1D-array of logical values with size(vec).
! If present, only those locations in vec corresponding to the true values in mask are used.
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! real(sp/dp) :: before (n-1)-th smallest value in input array, e.g. for median/percentile calculations
! real(sp/dp) :: previous same as before
! real(sp/dp) :: after (n+1)-th smallest value in input array
! real(sp/dp) :: next same as after
! RESTRICTIONS
! None
! EXAMPLE
! vec = (/ 1.,2.,3.,4.,5.,6.,7.,8.,9.,10. /)
! ! Returns 4
! out = n_element(vec,4)
! -> see also example in test directory
! LITERATURE
! Niklaus Wirth. "Algorithms and Data Structures". Prentice-Hall, Inc., 1985. ISBN 0-13-022005-1.
! HISTORY
! Written, Matthias Cuntz, May 2014 - based on qmedian
INTERFACE n_element
MODULE PROCEDURE n_element_dp, n_element_sp
END INTERFACE n_element
! ------------------------------------------------------------------
! NAME
! percentile
! PURPOSE
! Returns the value below which a certain percent of array values fall.
!
! If an optinal mask is given, values only on those locations that correspond
! to true values in the mask are used.
!
! Different definitions can be applied to interpolate the stepwise CDF of the given data.
! (1) Inverse empirical CDF (no interpolation, default MATHEMATICA)
! (2) Linear interpolation (California method)
! (3) Element numbered closest
! (4) Linear interpolation (hydrologist method)
! (5) Mean-based estimate (Weibull method, default IMSL)
! (6) Mode-based estimate
! (7) Median-based estimate
! (8) normal distribution estimate
!
! See: http://reference.wolfram.com/mathematica/tutorial/BasicStatistics.html
! CALLING SEQUENCE
! out = percentile(vec,k,mask=mask,mode_in=mode)
! INTENT(IN)
! real(sp/dp) :: vec(:) 1D-array with input numbers
! real(sp/dp) :: k[(:)] Percentage of percentile, can be 1 dimensional
! INTENT(INOUT)
! None
! INTENT(OUT)
! real(sp/dp) :: out[(size(k))] k-th percentile of values in input array, can be
! 1 dimensional corresponding to k
! INTENT(IN), OPTIONAL
! logical :: mask(:) 1D-array of logical values with size(vec).
! If present, only those locations in vec corresponding to the true values in mask are used.
! integer(i4) :: mode_in Specifies the interpolation scheme applied.
! Default:
! Inverse empirical CDF (no interpolation, default Mathematica)
! mode_in = 1_i4
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! vec = (/ 1.,2.,3.,4.,5.,6.,7.,8.,9.,10. /)
! ! Returns 10.
! out = percentile(vec,95.)
! ! Returns (10.,8)
! out = percentile(vec,(/95.,80./))
! -> see also example in test directory
! LITERATURE
! None
! HISTORY
! Written, Matthias Cuntz, Mar 2011
! Modified, Stephan Thober, Dec 2011 - added 1 dimensional version
! Juliane Mai, July 2012 - different interpolation schemes
! Modified, Matthias Cuntz, Juliane Mai, Jul 2012 - uses previous of ksmallest to half execution time
INTERFACE percentile
MODULE PROCEDURE percentile_0d_sp, percentile_0d_dp, percentile_1d_sp, percentile_1d_dp
END INTERFACE percentile
! ------------------------------------------------------------------
! NAME
! qmedian
! PURPOSE
! Quick calculation of the median thereby rearranging the input array.
! CALLING SEQUENCE
! out = qmedian(vec)
! INTENT(IN)
! None
! INTENT(INOUT)
! real(sp/dp) :: vec(:) 1D-array with input numbers.
! Wil be rearranged on output.
! INTENT(OUT)
! real(sp/dp) :: out median of values in input array
! INTENT(IN), OPTIONAL
! None
! INTENT(INOUT), OPTIONAL
! None
! INTENT(OUT), OPTIONAL
! None
! RESTRICTIONS
! None
! EXAMPLE
! vec = (/ 1.,2.,3.,4.,5.,6.,7.,8.,9.,10. /)
! ! Returns 5.5
! out = qmedian(vec)
! -> see also example in test directory
! LITERATURE
! Niklaus Wirth. "Algorithms and Data Structures". Prentice-Hall, Inc., 1985. ISBN 0-13-022005-1.
! HISTORY
! Written, Filip Hroch as part of Munipack: http://munipack.physics.muni.cz
! Modified, Matthias Cuntz, Jul 2012 - function, k=n/2+1
! Modified, Matthias Cuntz, Juliane Mai, Jul 2012 - real median for even n
INTERFACE qmedian
MODULE PROCEDURE qmedian_sp, qmedian_dp
END INTERFACE qmedian
! ------------------------------------------------------------------
PRIVATE
! ------------------------------------------------------------------
CONTAINS
! ------------------------------------------------------------------
FUNCTION median_dp(arrin, mask)
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: arrin
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
REAL(dp) :: median_dp
INTEGER(i4) :: n
REAL(dp), DIMENSION(:), ALLOCATABLE :: arr
REAL(dp) :: tmp
if (present(mask)) then
n = count(mask)
allocate(arr(n))
arr = pack(arrin,mask)
if (n < 2) stop 'median_dp: n < 2'
if (mod(n,2) == 0) then ! Even
median_dp = n_element(arr,n/2+1,previous=tmp)
median_dp = 0.5_dp*(median_dp+tmp)
else ! Odd
median_dp = n_element(arr,(n+1)/2)
endif
deallocate(arr)
else
n = size(arrin)
if (n < 2) stop 'median_dp: n < 2'
if (mod(n,2) == 0) then ! Even
median_dp = n_element(arrin,n/2+1,previous=tmp)
median_dp = 0.5_dp*(median_dp+tmp)
else ! Odd
median_dp = n_element(arrin,(n+1)/2)
endif
endif
END FUNCTION median_dp
FUNCTION median_sp(arrin, mask)
IMPLICIT NONE
REAL(sp), DIMENSION(:), INTENT(IN) :: arrin
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
REAL(sp) :: median_sp
INTEGER(i4) :: n
REAL(sp), DIMENSION(:), ALLOCATABLE :: arr
REAL(sp) :: tmp
if (present(mask)) then
n = count(mask)
allocate(arr(n))
arr = pack(arrin,mask)
if (n < 2) stop 'median_sp: n < 2'
if (mod(n,2) == 0) then ! Even
median_sp = n_element(arr,n/2+1,previous=tmp)
median_sp = 0.5_sp*(median_sp+tmp)
else ! Odd
median_sp = n_element(arr,(n+1)/2)
endif
deallocate(arr)
else
n = size(arrin)
if (n < 2) stop 'median_sp: n < 2'
if (mod(n,2) == 0) then ! Even
median_sp = n_element(arrin,n/2+1,previous=tmp)
median_sp = 0.5_sp*(median_sp+tmp)
else ! Odd
median_sp = n_element(arrin,(n+1)/2)
endif
endif
END FUNCTION median_sp
! ------------------------------------------------------------------
function n_element_dp(idat, n, mask, before, after, previous, next)
IMPLICIT NONE
real(dp), dimension(:), intent(in) :: idat
integer(i4), intent(in) :: n
logical, dimension(:), optional, intent(in) :: mask
real(dp), optional, intent(out) :: before
real(dp), optional, intent(out) :: after
real(dp), optional, intent(out) :: previous
real(dp), optional, intent(out) :: next
real(dp) :: n_element_dp
real(dp), dimension(:), allocatable :: dat
real(dp) :: w
integer(i4) :: nn, k
integer(i4) :: l,r,i,j
if (present(mask)) then
nn = count(mask)
allocate(dat(nn))
dat = pack(idat,mask)
else
nn = size(idat)
allocate(dat(nn))
dat = idat
endif
if (n < 1) stop 'n_element_dp: n < 1'
if (n > nn) stop 'n_element_dp: n > size(pack(dat,mask))'
dat = idat
nn = size(dat)
k = n !nn/2 + 1
l = 1
r = nn
do while( l < r )
n_element_dp = dat(k)
i = l
j = r
do
do while( dat(i) < n_element_dp )
i = i + 1
enddo
do while( n_element_dp < dat(j) )
j = j - 1
enddo
if ( i <= j ) then
w = dat(i)
dat(i) = dat(j)
dat(j) = w
i = i + 1
j = j - 1
endif
if ( i > j ) exit
enddo
if ( j < k ) l = i
if ( k < i ) r = j
enddo
! if (mod(nn,2) == 0) then
! n_element_dp = 0.5_dp*(dat(k) + maxval(dat(:k-1)))
! else
! n_element_dp = dat(k)
! endif
n_element_dp = dat(k)
if (present(before)) before = maxval(dat(:k-1))
if (present(previous)) previous = maxval(dat(:k-1))
if (present(after)) after = minval(dat(k+1:))
if (present(next)) next = minval(dat(k+1:))
deallocate(dat)
end function n_element_dp
function n_element_sp(idat, n, mask, before, after, previous, next)
IMPLICIT NONE
real(sp), dimension(:), intent(in) :: idat
integer(i4), intent(in) :: n
logical, dimension(:), optional, intent(in) :: mask
real(sp), optional, intent(out) :: before
real(sp), optional, intent(out) :: after
real(sp), optional, intent(out) :: previous
real(sp), optional, intent(out) :: next
real(sp) :: n_element_sp
real(sp), dimension(:), allocatable :: dat
real(sp) :: w
integer(i4) :: nn, k
integer(i4) :: l,r,i,j
if (present(mask)) then
nn = count(mask)
allocate(dat(nn))
dat = pack(idat,mask)
else
nn = size(idat)
allocate(dat(nn))
dat = idat
endif
if (n < 1) stop 'n_element_sp: n < 1'
if (n > nn) stop 'n_element_sp: n > size(pack(dat,mask))'
dat = idat
nn = size(dat)
k = n !nn/2 + 1
l = 1
r = nn
do while( l < r )
n_element_sp = dat(k)
i = l
j = r
do
do while( dat(i) < n_element_sp )
i = i + 1
enddo
do while( n_element_sp < dat(j) )
j = j - 1
enddo
if ( i <= j ) then
w = dat(i)
dat(i) = dat(j)
dat(j) = w
i = i + 1
j = j - 1
endif
if ( i > j ) exit
enddo
if ( j < k ) l = i
if ( k < i ) r = j
enddo
! if (mod(nn,2) == 0) then
! n_element_sp = 0.5_sp*(dat(k) + maxval(dat(:k-1)))
! else
! n_element_sp = dat(k)
! endif
n_element_sp = dat(k)
if (present(before)) before = maxval(dat(:k-1))
if (present(previous)) previous = maxval(dat(:k-1))
if (present(after)) after = minval(dat(k+1:))
if (present(next)) next = minval(dat(k+1:))
deallocate(dat)
end function n_element_sp
! ------------------------------------------------------------------
FUNCTION percentile_0d_dp(arrin,k,mask,mode_in)
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: arrin
REAL(dp), INTENT(IN) :: k
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
INTEGER(i4), OPTIONAL, INTENT(IN) :: mode_in
REAL(dp) :: percentile_0d_dp
INTEGER(i4) :: n, nn1, nn2
INTEGER(i4) :: mode
REAL(dp) :: kk, ks1, ks2
REAL(dp), DIMENSION(:), ALLOCATABLE :: arr
if (present(mask)) then
n = count(mask)
else
n = size(arrin)
endif
if (present(mode_in)) then
mode = mode_in
else
! Default : Inverse empirical CDF
mode = 1_i4
end if
if (n < 2) stop 'percentile_0d_dp: n < 2'
select case (mode)
! Inverse empirical CDF: Mathematica default
case(1_i4)
kk = k/100._dp*real(n,dp)
nn1 = min(n, max(1_i4,ceiling(kk,kind=i4)))
nn2 = nn1
! Linear interpolation (California method)
case(2_i4)
kk = k/100._dp*real(n,dp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Element numbered closest
case(3_i4)
kk = 0.5_dp+k/100._dp*real(n,dp)
nn1 = min(n, max(1_i4,floor(kk,kind=i4)))
nn2 = nn1
! Linear interpolation (hydrologist method)
case(4_i4)
kk = 0.5_dp+k/100._dp*(real(n,dp))
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Mean-based estimate (Weibull method): IMSL default
case(5_i4)
kk = k/100._dp*(real(n,dp)+1._dp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Mode-based estimate
case(6_i4)
kk = 1.0_dp+k/100._dp*(real(n,dp)-1._dp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Median-based estimate
case(7_i4)
kk = 1.0_dp/3.0_dp+k/100._dp*(real(n,dp)+1.0_dp/3.0_dp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Normal distribution estimate
case(8_i4)
kk = 3.0_dp/8.0_dp+k/100._dp*(real(n,dp)+1.0_dp/4.0_dp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! No valid mode
case default
stop 'percentile_0d_dp: mode > 8 not implemented'
end select
if (present(mask)) then
allocate(arr(n))
arr = pack(arrin,mask)
if (nn1 .eq. nn2) then
! no interpolation
percentile_0d_dp = n_element(arr,nn1)
else
! interpolation
ks2 = n_element(arr,nn2,previous=ks1)
percentile_0d_dp = ks1 + (ks2-ks1)*(kk-real(nn1,dp))
end if
deallocate(arr)
else
if (nn1 .eq. nn2) then
! no interpolation
percentile_0d_dp = n_element(arrin,nn1)
else
! interpolation
ks2 = n_element(arrin,nn2,previous=ks1)
percentile_0d_dp = ks1 + (ks2-ks1)*(kk-real(nn1,dp))
end if
endif
END FUNCTION percentile_0d_dp
FUNCTION percentile_0d_sp(arrin,k,mask,mode_in)
IMPLICIT NONE
REAL(sp), DIMENSION(:), INTENT(IN) :: arrin
REAL(sp), INTENT(IN) :: k
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
INTEGER(i4), OPTIONAL, INTENT(IN) :: mode_in
REAL(sp) :: percentile_0d_sp
INTEGER(i4) :: n, nn1, nn2
INTEGER(i4) :: mode
REAL(sp) :: kk, ks1, ks2
REAL(sp), DIMENSION(:), ALLOCATABLE :: arr
if (present(mask)) then
n = count(mask)
else
n = size(arrin)
endif
if (present(mode_in)) then
mode = mode_in
else
! Default : Inverse empirical CDF
mode = 1_i4
end if
if (n < 2) stop 'percentile_0d_sp: n < 2'
select case (mode)
! Inverse empirical CDF: Mathematica default
case(1_i4)
kk = k/100._sp*real(n,sp)
nn1 = min(n, max(1_i4,ceiling(kk,kind=i4)))
nn2 = nn1
! Linear interpolation (California method)
case(2_i4)
kk = k/100._sp*real(n,sp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Element numbered closest
case(3_i4)
kk = 0.5_sp+k/100._sp*real(n,sp)
nn1 = min(n, max(1_i4,floor(kk,kind=i4)))
nn2 = nn1
! Linear interpolation (hydrologist method)
case(4_i4)
kk = 0.5_sp+k/100._sp*(real(n,sp))
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Mean-based estimate (Weibull method): IMSL default
case(5_i4)
kk = k/100._sp*(real(n,sp)+1._sp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Mode-based estimate
case(6_i4)
kk = 1.0_sp+k/100._sp*(real(n,sp)-1._sp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Median-based estimate
case(7_i4)
kk = 1.0_sp/3.0_sp+k/100._sp*(real(n,sp)+1.0_sp/3.0_sp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! Normal distribution estimate
case(8_i4)
kk = 3.0_sp/8.0_sp+k/100._sp*(real(n,sp)+1.0_sp/4.0_sp)
nn1 = min(n, max(1_i4, floor(kk,kind=i4)))
nn2 = min(n, max(1_i4,ceiling(kk,kind=i4)))
! No valid mode
case default
stop 'percentile_0d_sp: mode > 8 not implemented'
end select
if (present(mask)) then
allocate(arr(n))
arr = pack(arrin,mask)
if (nn1 .eq. nn2) then
! no interpolation
percentile_0d_sp = n_element(arr,nn1)
else
! interpolation
ks2 = n_element(arr,nn2,previous=ks1)
percentile_0d_sp = ks1 + (ks2-ks1)*(kk-real(nn1,sp))
end if
deallocate(arr)
else
if (nn1 .eq. nn2) then
! no interpolation
percentile_0d_sp = n_element(arrin,nn1)
else
! interpolation
ks2 = n_element(arrin,nn2,previous=ks1)
percentile_0d_sp = ks1 + (ks2-ks1)*(kk-real(nn1,sp))
end if
endif
END FUNCTION percentile_0d_sp
function percentile_1d_dp(arrin,k,mask,mode_in)
IMPLICIT NONE
REAL(dp), DIMENSION(:), INTENT(IN) :: arrin
REAL(dp), DIMENSION(:), INTENT(IN) :: k
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
INTEGER(i4), OPTIONAL, INTENT(IN) :: mode_in
REAL(dp), DIMENSION(size(k)) :: percentile_1d_dp
INTEGER(i4) :: i, n
INTEGER(i4) :: mode
INTEGER(i4), DIMENSION(size(k)) :: nn1, nn2
REAL(dp), DIMENSION(size(k)) :: kk
REAL(dp) :: ks1, ks2
REAL(dp), DIMENSION(:), ALLOCATABLE :: arr
if (present(mask)) then
n = count(mask)
else
n = size(arrin)
endif
if (present(mode_in)) then
mode = mode_in
else
! Default : Inverse empirical CDF
mode = 1_i4
end if
! check consistency
!if (size(k) > size(arr)) stop 'percentile_1d_dp: more Quantiles than data: size(k) > size(arr)'
if (n < 2) stop 'percentile_1d_dp: n < 2'
select case (mode)
! Inverse empirical CDF: Mathematica default
case(1_i4)
kk(:) = k(:)/100._dp*real(n,dp)
nn1(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
nn2 = nn1
! Linear interpolation (California method)
case(2_i4)
kk(:) = k(:)/100._dp*real(n,dp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Element numbered closest
case(3_i4)
kk(:) = 0.5_dp+k(:)/100._dp*real(n,dp)
nn1(:) = min(n, max(1_i4,floor(kk(:),kind=i4)))
nn2 = nn1
! Linear interpolation (hydrologist method)
case(4_i4)
kk(:) = 0.5_dp+k(:)/100._dp*(real(n,dp))
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Mean-based estimate (Weibull method): IMSL default
case(5_i4)
kk(:) = k(:)/100._dp*(real(n,dp)+1._dp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Mode-based estimate
case(6_i4)
kk(:) = 1.0_dp+k(:)/100._dp*(real(n,dp)-1._dp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Median-based estimate
case(7_i4)
kk(:) = 1.0_dp/3.0_dp+k(:)/100._dp*(real(n,dp)+1.0_dp/3.0_dp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Normal distribution estimate
case(8_i4)
kk(:) = 3.0_dp/8.0_dp+k(:)/100._dp*(real(n,dp)+1.0_dp/4.0_dp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! No valid mode
case default
stop 'percentile_1d_dp: mode > 8 not implemented'
end select
if (present(mask)) then
allocate(arr(n))
arr = pack(arrin,mask)
do i=1, size(k)
if (nn1(i) .eq. nn2(i)) then
! no interpolation
percentile_1d_dp(i) = n_element(arr,nn1(i))
else
! interpolation
ks2 = n_element(arr,nn2(i),previous=ks1)
percentile_1d_dp(i) = ks1 + (ks2-ks1)*(kk(i)-real(nn1(i),dp))
end if
end do
deallocate(arr)
else
do i=1, size(k)
if (nn1(i) .eq. nn2(i)) then
! no interpolation
percentile_1d_dp(i) = n_element(arrin,nn1(i))
else
! interpolation
ks2 = n_element(arrin,nn2(i),previous=ks1)
percentile_1d_dp(i) = ks1 + (ks2-ks1)*(kk(i)-real(nn1(i),dp))
end if
end do
endif
END function percentile_1d_dp
function percentile_1d_sp(arrin,k,mask,mode_in)
IMPLICIT NONE
REAL(sp), DIMENSION(:), INTENT(IN) :: arrin
REAL(sp), DIMENSION(:), INTENT(IN) :: k
LOGICAL, DIMENSION(:), OPTIONAL, INTENT(IN) :: mask
INTEGER(i4), OPTIONAL, INTENT(IN) :: mode_in
REAL(sp), DIMENSION(size(k)) :: percentile_1d_sp
INTEGER(i4) :: i, n
INTEGER(i4) :: mode
INTEGER(i4), DIMENSION(size(k)) :: nn1, nn2
REAL(sp), DIMENSION(size(k)) :: kk
REAL(sp) :: ks1, ks2
REAL(sp), DIMENSION(:), ALLOCATABLE :: arr
if (present(mask)) then
n = count(mask)
else
n = size(arrin)
endif
if (present(mode_in)) then
mode = mode_in
else
! Default : Inverse empirical CDF
mode = 1_i4
end if
! check consistency
!if (size(k) > size(arr)) stop 'percentile_1d_sp: more Quantiles than data: size(k) > size(arr)'
if (n < 2) stop 'percentile_1d_sp: n < 2'
select case (mode)
! Inverse empirical CDF: Mathematica default
case(1_i4)
kk(:) = k(:)/100._sp*real(n,sp)
nn1(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
nn2 = nn1
! Linear interpolation (California method)
case(2_i4)
kk(:) = k(:)/100._sp*real(n,sp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Element numbered closest
case(3_i4)
kk(:) = 0.5_sp+k(:)/100._sp*real(n,sp)
nn1(:) = min(n, max(1_i4,floor(kk(:),kind=i4)))
nn2 = nn1
! Linear interpolation (hydrologist method)
case(4_i4)
kk(:) = 0.5_sp+k(:)/100._sp*(real(n,sp))
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Mean-based estimate (Weibull method): IMSL default
case(5_i4)
kk(:) = k(:)/100._sp*(real(n,sp)+1._sp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Mode-based estimate
case(6_i4)
kk(:) = 1.0_sp+k(:)/100._sp*(real(n,sp)-1._sp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Median-based estimate
case(7_i4)
kk(:) = 1.0_sp/3.0_sp+k(:)/100._sp*(real(n,sp)+1.0_sp/3.0_sp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! Normal distribution estimate
case(8_i4)
kk(:) = 3.0_sp/8.0_sp+k(:)/100._sp*(real(n,sp)+1.0_sp/4.0_sp)
nn1(:) = min(n, max(1_i4, floor(kk(:),kind=i4)))
nn2(:) = min(n, max(1_i4,ceiling(kk(:),kind=i4)))
! No valid mode
case default
stop 'percentile_1d_sp: mode > 8 not implemented'
end select
if (present(mask)) then
allocate(arr(n))
arr = pack(arrin,mask)
do i=1, size(k)
if (nn1(i) .eq. nn2(i)) then
! no interpolation
percentile_1d_sp(i) = n_element(arr,nn1(i))
else
! interpolation
ks2 = n_element(arr,nn2(i),previous=ks1)
percentile_1d_sp(i) = ks1 + (ks2-ks1)*(kk(i)-real(nn1(i),sp))
end if
end do
deallocate(arr)
else
do i=1, size(k)
if (nn1(i) .eq. nn2(i)) then
! no interpolation
percentile_1d_sp(i) = n_element(arrin,nn1(i))
else
! interpolation
ks2 = n_element(arrin,nn2(i),previous=ks1)
percentile_1d_sp(i) = ks1 + (ks2-ks1)*(kk(i)-real(nn1(i),sp))
end if
end do
endif
END function percentile_1d_sp
! ------------------------------------------------------------------
function qmedian_dp(dat)
IMPLICIT NONE
real(dp), dimension(:), intent(inout) :: dat
real(dp) :: qmedian_dp
real(dp) :: w
integer(i4) :: n,k
integer(i4) :: l,r,i,j
n = size(dat)
k = n/2 + 1
l = 1
r = n
do while( l < r )