-
Notifications
You must be signed in to change notification settings - Fork 375
/
Copy pathgmt_stat.c
2907 lines (2527 loc) · 91.8 KB
/
gmt_stat.c
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
/*--------------------------------------------------------------------
*
* Copyright (c) 1991-2025 by the GMT Team (https://www.generic-mapping-tools.org/team.html)
* See LICENSE.TXT file for copying and redistribution conditions.
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published by
* the Free Software Foundation; version 3 or any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details.
*
* Contact info: www.generic-mapping-tools.org
*--------------------------------------------------------------------*/
/*
* Misc statistical and special functions.
*
* Author: Walter H. F. Smith, P. Wessel, R. Scharroo
* Date: 1-JAN-2010
* Version: 5.x
*
* A) List of exported gmt_* functions available to modules and libraries via gmt_dev.h:
*
* gmt_Fcrit
* gmt_PvQv
* gmt_bei
* gmt_ber
* gmt_binom_cdf
* gmt_binom_pdf
* gmt_chebyshev
* gmt_chi2
* gmt_chi2_pdf
* gmt_chi2crit
* gmt_combination
* gmt_corrcoeff
* gmt_corrcoeff_f
* gmt_dilog
* gmt_erfinv
* gmt_extreme
* gmt_f_cdf
* gmt_f_pdf
* gmt_factorial
* gmt_fisher_pdf
* gmt_get_cellarea
* gmt_getmad
* gmt_getmad_f
* gmt_grd_lmsscl
* gmt_grd_mad
* gmt_grd_mean
* gmt_grd_median
* gmt_grd_mode
* gmt_grd_rms
* gmt_grd_std
* gmt_i0
* gmt_i1
* gmt_in
* gmt_k0
* gmt_k1
* gmt_kei
* gmt_ker
* gmt_kn
* gmt_lrand
* gmt_mean_and_std
* gmt_mean_weighted
* gmt_median
* gmt_median_weighted
* gmt_mode
* gmt_mode_f
* gmt_mode_weighted
* gmt_nrand
* gmt_permutation
* gmt_plm
* gmt_plm_bar
* gmt_plm_bar_all
* gmt_poisson_cdf
* gmt_poissonpdf
* gmt_psi
* gmt_quantile
* gmt_quantile_f
* gmt_quantile_weighted
* gmt_rand
* gmt_sig_f
* gmt_sinc
* gmt_std_weighted
* gmt_t_cdf
* gmt_t_pdf
* gmt_tcrit
* gmt_von_mises_mu_and_kappa
* gmt_vonmises_pdf
* gmt_weibull_cdf
* gmt_weibull_crit
* gmt_weibull_pdf
* gmt_zcrit
* gmt_zdist
*
* B) List of exported gmtlib_* functions available to libraries via gmt_internals.h:
*
* gmtlib_compare_observation
*/
#include "gmt_dev.h"
#include "gmt_internals.h"
enum gmtstat_enum_cplx {GMT_RE = 0, GMT_IM = 1}; /* Real and imaginary indices */
/* --------- Local functions to this file ------- */
#define ITMAX 100
GMT_LOCAL double gmtstat_ln_gamma (struct GMT_CTRL *GMT, double xx) {
/* Routine to compute natural log of Gamma(x)
by Lanczos approximation. Most accurate
for x > 1; fails for x <= 0. No error
checking is done here; it is assumed
that this is called by gmtstat_ln_gamma_r() */
static double cof[6] = {
76.18009173,
-86.50532033,
24.01409822,
-1.231739516,
0.120858003e-2,
-0.536382e-5
};
double x, tmp, ser;
int i;
x = xx - 1.0;
tmp = x + 5.5;
tmp = (x + 0.5) * d_log (GMT,tmp) - tmp;
ser = 1.0;
for (i = 0; i < 6; i++) {
x += 1.0;
ser += (cof[i]/x);
}
return (tmp + d_log (GMT, 2.50662827465*ser) );
}
GMT_LOCAL int gmtstat_ln_gamma_r (struct GMT_CTRL *GMT, double x, double *lngam) {
/* Get natural logarithm of Gamma(x), x > 0.
To maintain full accuracy, this
routine uses Gamma(1 + x) / x when
x < 1. This routine in turn calls
gmtstat_ln_gamma(x), which computes the
actual function value. gmtstat_ln_gamma
assumes it is being called in a
smart way, and does not check the
range of x. */
if (x > 1.0) {
*lngam = gmtstat_ln_gamma (GMT, x);
return (0);
}
if (x > 0.0 && x < 1.0) {
*lngam = gmtstat_ln_gamma (GMT, 1.0 + x) - d_log (GMT,x);
return (0);
}
if (x == 1.0) {
*lngam = 0.0;
return (0);
}
GMT_Report (GMT->parent, GMT_MSG_WARNING, "Ln Gamma: Bad x (x <= 0).\n");
return (GMT_NOTSET);
}
GMT_LOCAL void gmtstat_gamma_ser (struct GMT_CTRL *GMT, double *gamser, double a, double x, double *gln) {
/* Returns the incomplete gamma function P(a,x) by series rep.
* Press et al, gser() */
int n;
double sum, del, ap;
if (gmtstat_ln_gamma_r (GMT, a, gln) == GMT_NOTSET) {
*gamser = GMT->session.d_NaN;
return;
}
if (x < 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "x < 0 in gmtstat_gamma_ser(x)\n");
*gamser = GMT->session.d_NaN;
return;
}
if (x == 0.0) {
*gamser = 0.0;
return;
}
ap = a;
del = sum = 1.0 / a;
for (n = 1; n <= ITMAX; n++) {
ap += 1.0;
del *= x / ap;
sum += del;
if (fabs (del) < fabs (sum) * DBL_EPSILON) {
*gamser = sum * exp (-x + a * log (x) - (*gln));
return;
}
}
GMT_Report (GMT->parent, GMT_MSG_WARNING, "a too large, ITMAX too small in gmtstat_gamma_ser(x)\n");
}
GMT_LOCAL void gmtstat_gamma_cf (struct GMT_CTRL *GMT, double *gammcf, double a, double x, double *gln) {
/* Returns the incomplete gamma function P(a,x) by continued fraction.
* Press et al, gcf() */
int n;
double gold = 0.0, g, fac = 1.0, b1 = 1.0;
double b0 = 0.0, anf, ana, an, a1, a0 = 1.0;
if (gmtstat_ln_gamma_r (GMT, a, gln) == GMT_NOTSET) {
*gln = GMT->session.d_NaN;
return;
}
a1 = x;
for (n = 1; n <= ITMAX; n++) {
an = (double) n;
ana = an - a;
a0 = (a1 + a0 * ana) * fac;
b0 = (b1 + b0 * ana) * fac;
anf = an * fac;
a1 = x * a0 + anf * a1;
b1 = x * b0 + anf * b1;
if (a1 != 0.0) {
fac = 1.0 / a1;
g = b1 * fac;
if (fabs ((g - gold) / g) < DBL_EPSILON) {
*gammcf = exp (-x + a * log (x) - (*gln)) * g;
return;
}
gold = g;
}
}
GMT_Report (GMT->parent, GMT_MSG_WARNING, "a too large, ITMAX too small in gmtstat_gamma_cf(x)\n");
}
GMT_LOCAL double gmtstat_gammq (struct GMT_CTRL *GMT, double a, double x) {
/* Returns Q(a,x) = 1 - P(a,x) Inc. Gamma function */
double G = 0.0, gln;
if (x < 0.0 || a <= 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "Invalid arguments to GMT_gammaq\n");
return (GMT->session.d_NaN);
}
if (x < (a + 1.0)) {
gmtstat_gamma_ser (GMT, &G, a, x, &gln);
return (1.0 - G);
}
gmtstat_gamma_cf (GMT, &G, a, x, &gln);
return (G);
}
#define BETA_EPS 3.0e-7
GMT_LOCAL double gmtstat_cf_beta (struct GMT_CTRL *GMT, double a, double b, double x) {
/* Continued fraction method called by gmtstat_inc_beta. */
double am = 1.0, bm = 1.0, az = 1.0;
double qab, qap, qam, bz, em, tem, d;
double ap, bp, app, bpp, aold;
int m = 0;
qab = a + b;
qap = a + 1.0;
qam = a - 1.0;
bz = 1.0 - qab * x / qap;
do {
m++;
em = (double)m;
tem = em + em;
d = em*(b-m)*x/((qam+tem)*(a+tem));
ap = az+d*am;
bp = bz+d*bm;
d = -(a+m)*(qab+em)*x/((a+tem)*(qap+tem));
app = ap+d*az;
bpp = bp+d*bz;
aold = az;
am = ap/bpp;
bm = bp/bpp;
az = app/bpp;
bz = 1.0;
} while (((fabs (az-aold) ) >= (BETA_EPS * fabs (az))) && (m < ITMAX));
if (m == ITMAX) GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_cf_beta: A or B too big, or ITMAX too small.\n");
return (az);
}
GMT_LOCAL int gmtstat_inc_beta (struct GMT_CTRL *GMT, double a, double b, double x, double *ibeta) {
double bt, gama, gamb, gamab;
if (a <= 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_inc_beta: Bad a (a <= 0).\n");
return(GMT_NOTSET);
}
if (b <= 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_inc_beta: Bad b (b <= 0).\n");
return(GMT_NOTSET);
}
if (x > 0.0 && x < 1.0) {
gmtstat_ln_gamma_r(GMT, a, &gama);
gmtstat_ln_gamma_r(GMT, b, &gamb);
gmtstat_ln_gamma_r(GMT, (a+b), &gamab);
bt = exp(gamab - gama - gamb
+ a * d_log (GMT, x) + b * d_log (GMT, 1.0 - x) );
/* Here there is disagreement on the range of x which
converges efficiently. Abramowitz and Stegun
say to use x < (a - 1) / (a + b - 2). Editions
of Numerical Recipes through mid 1987 say
x < ( (a + 1) / (a + b + 1), but the code has
x < ( (a + 1) / (a + b + 2). Editions printed
late 1987 and after say x < ( (a + 1) / (a + b + 2)
in text as well as code. What to do ? */
if (x < ( (a + 1) / (a + b + 2) ) )
*ibeta = bt * gmtstat_cf_beta (GMT, a, b, x) / a;
else
*ibeta = 1.0 - bt * gmtstat_cf_beta (GMT, b, a, (1.0 - x) ) / b;
return(0);
}
else if (x == 0.0) {
*ibeta = 0.0;
return (0);
}
else if (x == 1.0) {
*ibeta = 1.0;
return (0);
}
else if (x < 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_inc_beta: Bad x (x < 0).\n");
*ibeta = 0.0;
}
else if (x > 1.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_inc_beta: Bad x (x > 1).\n");
*ibeta = 1.0;
}
return (GMT_NOTSET);
}
GMT_LOCAL int gmtstat_f_q (struct GMT_CTRL *GMT, double chisq1, uint64_t nu1, double chisq2, uint64_t nu2, double *prob) {
/* Routine to compute Q(F, nu1, nu2) = 1 - P(F, nu1, nu2), where nu1
and nu2 are positive integers, chisq1 and chisq2 are random
variables having chi-square distributions with nu1 and nu2
degrees of freedom, respectively (chisq1 and chisq2 >= 0.0),
F = (chisq1/nu1)/(chisq2/nu2) has the F-distribution, and
P(F, nu1, nu2) is the cumulative F-distribution, that is,
the integral from 0 to (chisq1/nu1)/(chisq2/nu2) of the F-
distribution. Q = 1 - P is small when (chisq1/nu1)/(chisq2/nu2)
is large with respect to 1. That is, the value returned by
this routine is the likelihood that an F >= (chisq1/nu1)/
(chisq2/nu2) would occur by chance.
Follows Abramowitz and Stegun.
This is different from the method in Numerical Recipes, which
uses the incomplete beta function but makes no use of the fact
that nu1 and nu2 are known to be integers, and thus there is
a finite limit on the sum for their expression.
W H F Smith, August, 1999.
REVISED by W H F Smith, October 27, 2000 after GMT 3.3.6 release.
I found that the A&S methods overflowed for large nu1 and nu2, so
I decided to go back to the gmtstat_inc_beta way of doing things.
*/
/* Check range of arguments: */
if (nu1 <= 0 || nu2 <= 0 || chisq1 < 0.0 || chisq2 < 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_f_q: Bad argument(s).\n");
return (GMT_NOTSET);
}
/* Extreme cases evaluate immediately: */
if (chisq1 == 0.0) {
*prob = 1.0;
return (0);
}
if (chisq2 == 0.0) {
*prob = 0.0;
return (0);
}
/* REVISION of Oct 27, 2000: This inc beta call here returns
the value. All subsequent code is not used. */
if (gmtstat_inc_beta (GMT, 0.5*nu2, 0.5*nu1, chisq2/(chisq2+chisq1), prob) ) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_f_q: Trouble in gmtstat_inc_beta call.\n");
return (GMT_NOTSET);
}
return (0);
}
GMT_LOCAL int gmtstat_f_test_new (struct GMT_CTRL *GMT, double chisq1, uint64_t nu1, double chisq2, uint64_t nu2, double *prob, int iside) {
/* Given chisq1 and chisq2, random variables distributed as chi-square
with nu1 and nu2 degrees of freedom, respectively, except that
chisq1 is scaled by var1, and chisq2 is scaled by var2, let
the null hypothesis, H0, be that var1 = var2. This routine
assigns prob, the probability that we can reject H0 in favor
of a new hypothesis, H1, according to iside:
iside=+1 means H1 is that var1 > var2
iside=-1 means H1 is that var1 < var2
iside=0 means H1 is that var1 != var2.
This routine differs from the old gmtstat_f_test() by adding the
argument iside and allowing one to choose the test. The old
routine in effect always set iside=0.
This routine also differs from gmtstat_f_test() in that the former
used the incomplete beta function and this one uses gmtstat_f_q().
Returns 0 on success, -1 on failure.
WHF Smith, 12 August 1999.
*/
double q; /* The probability from gmtstat_f_q(), which is the prob
that H0 should be retained even though
chisq1/nu1 > chisq2/nu2. */
if (chisq1 <= 0.0 || chisq2 <= 0.0 || nu1 < 1 || nu2 < 1) {
*prob = GMT->session.d_NaN;
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_f_test_new: Bad argument(s).\n");
return (GMT_NOTSET);
}
gmtstat_f_q (GMT, chisq1, nu1, chisq2, nu2, &q);
if (iside > 0)
*prob = 1.0 - q;
else if (iside < 0)
*prob = q;
else if ((chisq1/nu1) <= (chisq2/nu2))
*prob = 2.0*q;
else
*prob = 2.0*(1.0 - q);
return (0);
}
GMT_LOCAL double gmtstat_factln (struct GMT_CTRL *GMT, int n) {
/* returns log(n!) */
static double a[101]; /* Automatically initialized to zero */
if (n < 0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "n < 0 in gmtstat_factln(n)\n");
return (GMT->session.d_NaN);
}
if (n <= 1) return 0.0;
if (n <= 100) return (a[n] ? a[n] : (a[n] = gmtstat_ln_gamma (GMT, n+1.0)));
else return gmtstat_ln_gamma (GMT, n+1.0);
}
GMT_LOCAL double gmtstat_beta (struct GMT_CTRL *GMT, double z, double w) {
double g1 = 0.0, g2 = 0.0, g3 = 0.0;
gmtstat_ln_gamma_r (GMT, z, &g1);
gmtstat_ln_gamma_r (GMT, w, &g2);
gmtstat_ln_gamma_r (GMT, z+w, &g3);
return exp (g1 + g2 - g3);
}
#if 0 /* Legacy code */
GMT_LOCAL int gmtstat_f_test (struct GMT_CTRL *GMT, double chisq1, uint64_t nu1, double chisq2, uint64_t nu2, double *prob) {
/* Routine to compute the probability that
two variances are the same.
chisq1 is distributed as chisq with
nu1 degrees of freedom; ditto for
chisq2 and nu2. If these are independent
and we form the ratio
F = max(chisq1,chisq2)/min(chisq1,chisq2)
then we can ask what is the probability
that an F greater than this would occur
by chance. It is this probability that
is returned in prob. When prob is small,
it is likely that the two chisq represent
two different populations; the confidence
that the two do not represent the same pop
is 1.0 - prob. This is a two-sided test.
This follows some ideas in Numerical Recipes, CRC Handbook,
and Abramowitz and Stegun. */
double f, df1, df2, p1, p2;
if (chisq1 <= 0.0) {
GMT_Report (GMT->parent, GMT_MSG_ERROR, "gmtstat_f_test: Chi-Square One <= 0.0\n");
return(GMT_NOTSET);
}
if (chisq2 <= 0.0) {
GMT_Report (GMT->parent, GMT_MSG_ERROR, "gmtstat_f_test: Chi-Square Two <= 0.0\n");
return(GMT_NOTSET);
}
if (chisq1 > chisq2) {
f = chisq1/chisq2;
df1 = (double)nu1;
df2 = (double)nu2;
}
else {
f = chisq2/chisq1;
df1 = (double)nu2;
df2 = (double)nu1;
}
if (gmtstat_inc_beta(GMT, 0.5*df2, 0.5*df1, df2/(df2+df1*f), &p1) ) {
GMT_Report (GMT->parent, GMT_MSG_ERROR, "gmtstat_f_test: Trouble on 1st gmtstat_inc_beta call.\n");
return(GMT_NOTSET);
}
if (gmtstat_inc_beta(GMT, 0.5*df1, 0.5*df2, df1/(df1+df2/f), &p2) ) {
GMT_Report (GMT->parent, GMT_MSG_ERROR, "gmtstat_f_test: Trouble on 2nd gmtstat_inc_beta call.\n");
return(GMT_NOTSET);
}
*prob = p1 + (1.0 - p2);
return (0);
}
#endif
GMT_LOCAL int gmtstat_student_t_a (struct GMT_CTRL *GMT, double t, uint64_t n, double *prob) {
/* Probability integral called A(t,n) by Abramowitz &
Stegun for the student's t distribution with n degrees
of freedom. Uses expressions A&S 26.7.3 and 26.7.4
If X is distributed N(0,1) and V is distributed chi-
square with n degrees of freedom, then
tau = X / sqrt(V/n) is said to have Student's t-
distribution with n degrees of freedom. For example,
tau could be the sample mean divided by the sample
standard deviation, for a sample of N points; then
n = N - 1.
This function sets *prob = GMT->session.d_NaN and returns (-1)
if t < 0. Otherwise it sets *prob = the probability
fabs(tau) <= t and returns (0).
As n -> oo, we can replace this function with
erf (t / M_SQRT2). However, it isn't clear how large
n has to be to make this a good approximation. I
consulted six books; one of them suggested this
approximation for n >= 30, but all the others did not
say when to use this approximation (A&S, in particular,
does not say). I tried some numerical experiments
which suggested that the relative error in this
approximation would be < 0.01 for n > 30, all t, but
I also found that the expression here is stable to
large n and large t, so I decided to leave it as is.
W H F Smith, August 1999.
*/
double theta, s, c, csq, term, sum;
int64_t k, kstop, kt, kb;
if (t < 0.0 || n == 0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmtstat_student_t_a: Bad argument(s).\n");
*prob = GMT->session.d_NaN;
return (GMT_NOTSET);
}
if (t == 0.0) {
*prob = 0.0;
return (0);
}
theta = atan (t/sqrt ((double)n));
if (n == 1) {
*prob = 2.0 * theta / M_PI;
return (0);
}
sincos (theta, &s, &c);
csq = c * c;
kstop = n-2;
if (n%2 == 1) {
kt = 0;
kb = 1;
k = 1;
term = c;
}
else {
kt = -1;
kb = 0;
k = 0;
term = 1.0;
}
sum = term;
while (k < kstop) {
k += 2;
kt += 2;
kb += 2;
term *= (kt * csq)/kb;
sum += term;
}
sum *= s;
if (n%2 == 1)
*prob = 2.0 * (theta + sum) / M_PI;
else
*prob = sum;
/* Adjust in case of roundoff: */
if (*prob < 0.0) *prob = 0.0;
if (*prob > 1.0) *prob = 1.0;
return (0);
}
GMT_LOCAL void gmtstat_Cmul (double A[], double B[], double C[]) {
/* Complex multiplication */
C[GMT_RE] = A[GMT_RE]*B[GMT_RE] - A[GMT_IM]*B[GMT_IM];
C[GMT_IM] = A[GMT_RE]*B[GMT_IM] + A[GMT_IM]*B[GMT_RE];
}
GMT_LOCAL void gmtstat_Cdiv (double A[], double B[], double C[]) {
/* Complex division */
double denom;
denom = B[GMT_RE]*B[GMT_RE] + B[GMT_IM]*B[GMT_IM];
C[GMT_RE] = (A[GMT_RE]*B[GMT_RE] + A[GMT_IM]*B[GMT_IM])/denom;
C[GMT_IM] = (A[GMT_IM]*B[GMT_RE] - A[GMT_RE]*B[GMT_IM])/denom;
}
#if 0 /* Unused */
GMT_LOCAL void gmtstat_Ccot (double Z[], double cotZ[]) {
/* Complex cot(z) */
double sx, cx, e, A[2], B[2];
sincos (2.0*Z[0], &sx, &cx);
e = exp (-2.0*Z[1]);
A[0] = -e * sx; A[1] = B[0] = e * cx;
A[1] += 1.0; B[0] -= 1.0; B[1] = -A[0];
gmtstat_Cdiv (A, B, cotZ);
}
#endif
GMT_LOCAL double gmtstat_Cabs (double A[]) {
return (hypot (A[GMT_RE], A[GMT_IM]));
}
GMT_LOCAL void gmtstat_F_to_ch1_ch2 (struct GMT_CTRL *GMT, double F, double nu1, double nu2, double *chisq1, double *chisq2) {
/* Silly routine to break F up into parts needed for gmtstat_f_q */
gmt_M_unused(GMT);
*chisq2 = 1.0;
*chisq1 = F * nu1 / nu2;
}
int gmtlib_compare_observation (const void *a, const void *b) {
const struct GMT_OBSERVATION *obs_1 = a, *obs_2 = b;
/* Sorts observations into ascending order based on obs->value */
if (obs_1->value < obs_2->value)
return -1;
if (obs_1->value > obs_2->value)
return 1;
return 0;
}
int gmt_sig_f (struct GMT_CTRL *GMT, double chi1, uint64_t n1, double chi2, uint64_t n2, double level, double *prob) {
/* Returns true if chi1/n1 significantly less than chi2/n2
at the level level. Returns false if:
error occurs in gmtstat_f_test_new();
chi1/n1 not significantly < chi2/n2 at level.
Changed 12 August 1999 to use gmtstat_f_test_new() */
int trouble;
trouble = gmtstat_f_test_new (GMT, chi1, n1, chi2, n2, prob, -1);
if (trouble) return (0);
return ((*prob) >= level);
}
/*
* Kelvin-Bessel functions, ber, bei, ker, kei. ber(x) and bei(x) are even;
* ker(x) and kei(x) are defined only for x > 0 and x >=0, respectively.
* For x <= 8 we use polynomial approximations in Abramowitz & Stegun.
* For x > 8 we use asymptotic series of Russell, quoted in Watson (Theory
* of Bessel Functions).
*/
double gmt_ber (struct GMT_CTRL *GMT, double x) {
double t, rxsq, alpha, beta;
gmt_M_unused(GMT);
if (x == 0.0) return (1.0);
/* ber is an even function of x: */
x = fabs (x);
if (x <= 8.0) {
/* Telescoped power series from Abramowitz & Stegun */
t = x * 0.125;
t *= t;
t *= t; /* t = pow(x/8, 4) */
return (1.0 + t*(-64.0 + t*(113.77777774 + t*(-32.36345652 + t*(2.64191397 + t*(-0.08349609 + t*(0.00122552 - 0.00000901 * t)))))));
}
else {
/* Russell's asymptotic approximation, from Watson, p. 204 */
rxsq = 1.0 / (x * x);
t = x / M_SQRT2;
alpha = t;
beta = t - 0.125 * M_PI;
t *= 0.125 * rxsq;
alpha += t;
beta -= t;
beta -= 0.0625*rxsq;
t *= (25.0/48.0)*rxsq;
alpha -= t;
beta -= t;
alpha -= (13.0/128.0)*(rxsq*rxsq);
return (exp (alpha) * cos (beta) / sqrt (2.0 * M_PI * x) );
}
}
double gmt_bei (struct GMT_CTRL *GMT, double x) {
double t, rxsq, alpha, beta;
gmt_M_unused(GMT);
if (x == 0.0) return (0.0);
/* bei is an even function of x: */
x = fabs(x);
if (x <= 8.0) {
/* Telescoped power series from Abramowitz & Stegun */
t = x * 0.125;
rxsq = t*t;
t = rxsq * rxsq; /* t = pow(x/8, 4) */
return (rxsq * (16.0 + t*(-113.77777774 + t*(72.81777742 + t*(-10.56765779 + t*(0.52185615 + t*(-0.01103667 + t*(0.00011346))))))));
}
else {
/* Russell's asymptotic approximation, from Watson, p. 204 */
rxsq = 1.0 / (x * x);
t = x / M_SQRT2;
alpha = t;
beta = t - 0.125 * M_PI;
t *= 0.125 * rxsq;
alpha += t;
beta -= t;
beta -= 0.0625*rxsq;
t *= (25.0/48.0)*rxsq;
alpha -= t;
beta -= t;
alpha -= (13.0/128.0)*(rxsq*rxsq);
return (exp (alpha) * sin (beta) / sqrt (2.0 * M_PI * x));
}
}
double gmt_ker (struct GMT_CTRL *GMT, double x) {
double t, rxsq, alpha, beta;
if (x <= 0.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "x <= 0 in gmt_ker(x)\n");
return (GMT->session.d_NaN);
}
if (x <= 8.0) {
/* Telescoped power series from Abramowitz & Stegun */
t = 0.125 * x;
t *= t;
t *= t; /* t = pow(x/8, 4) */
return (-log (0.5 * x) * gmt_ber (GMT, x) + 0.25 * M_PI * gmt_bei (GMT, x) -M_EULER + \
t * (-59.05819744 + t * (171.36272133 + t * (-60.60977451 + t * (5.65539121 + t * (-0.199636347 + t * (0.00309699 + t * (-0.00002458 * t))))))));
}
else {
/* Russell's asymptotic approximation, from Watson, p. 204 */
rxsq = 1.0 / (x * x);
t = -x / M_SQRT2;
alpha = t;
beta = t - 0.125 * M_PI;
t *= 0.125 * rxsq;
alpha += t;
beta -= t;
beta -= 0.0625*rxsq;
t *= (25.0/48.0)*rxsq;
alpha -= t;
beta -= t;
alpha -= (13.0/128.0)*(rxsq*rxsq);
return (exp (alpha) * cos (beta) / sqrt (2.0 * x / M_PI));
}
}
double gmt_kei (struct GMT_CTRL *GMT, double x) {
double t, rxsq, alpha, beta;
if (x <= 0.0) {
/* Zero is valid. If near enough to zero, return kei(0) */
if (x > -GMT_CONV8_LIMIT) return (-0.25 * M_PI);
GMT_Report (GMT->parent, GMT_MSG_WARNING, "x < 0 in gmt_kei(x)\n");
return (GMT->session.d_NaN);
}
if (x <= 8.0) {
/* Telescoped power series from Abramowitz & Stegun */
t = x * 0.125;
rxsq = t*t;
t = rxsq * rxsq; /* t = pow(x/8, 4) */
return (-log (0.5 * x) * gmt_bei (GMT, x) - 0.25 * M_PI * gmt_ber (GMT, x) +
rxsq * (6.76454936 + t * (-142.91827687 + t * (124.23569650 + t * (-21.30060904 + t * (1.17509064 + t * (-0.02695875 + t * (0.00029532 * t))))))));
}
else {
/* Russell's asymptotic approximation, from Watson, p. 204 */
rxsq = 1.0 / (x * x);
t = -x / M_SQRT2;
alpha = t;
beta = t - 0.125 * M_PI;
t *= 0.125 * rxsq;
alpha += t;
beta -= t;
beta -= 0.0625*rxsq;
t *= (25.0/48.0)*rxsq;
alpha -= t;
beta -= t;
alpha -= (13.0/128.0)*(rxsq*rxsq);
return (exp (alpha) * sin (beta) / sqrt (2.0 * x / M_PI));
}
}
double gmt_i0 (struct GMT_CTRL *GMT, double x) {
/* Modified from code in Press et al. */
double y, res;
gmt_M_unused(GMT);
if (x < 0.0) x = -x;
if (x < 3.75) {
y = x * x / 14.0625;
res = 1.0 + y * (3.5156229 + y * (3.0899424 + y * (1.2067492 + y * (0.2659732 + y * (0.360768e-1 + y * 0.45813e-2)))));
}
else {
y = 3.75 / x;
res = (exp (x) / sqrt (x)) * (0.39894228 + y * (0.1328592e-1 + y * (0.225319e-2 + y * (-0.157565e-2 + y * (0.916281e-2 + y * (-0.2057706e-1 + y * (0.2635537e-1 + y * (-0.1647633e-1 + y * 0.392377e-2))))))));
}
return (res);
}
double gmt_i1 (struct GMT_CTRL *GMT, double x) {
/* Modified Bessel function I1(x) */
double y, res;
gmt_M_unused(GMT);
if (x < 0.0) x = -x;
if (x < 3.75) {
y = pow (x / 3.75, 2.0);
res = x * (0.5 + y * (0.87890594 + y * (0.51498869 + y * (0.15084934 + y * (0.02658733 + y * (0.00301532 + y * 0.00032411))))));
}
else {
y = 3.75 / x;
res = (exp (x) / sqrt (x)) * (0.39894228 + y * (-0.03988024 + y * (-0.00362018 + y * (0.00163801+ y * (-0.01031555 + y * (0.02282967 + y * (-0.02895312 + y * (0.01787654 -y * 0.00420059))))))));
if (x < 0.0) res = - res;
}
return (res);
}
double gmt_in (struct GMT_CTRL *GMT, unsigned int n, double x) {
/* Modified Bessel function In(x) */
unsigned int j, m, IACC = 40;
double res, tox, bip, bi, bim;
double BIGNO = 1.0e10, BIGNI = 1.0e-10;
if (n == 0) return (gmt_i0 (GMT, x));
if (n == 1) return (gmt_i1 (GMT, x));
if (x == 0.0) return (0.0);
tox = 2.0 / fabs (x);
bip = res = 0.0;
bi = 1.0;
m = 2 * (n + urint (sqrt ((double)(IACC * n))));
for (j = m; j >= 1; j--) {
bim = bip + ((double)j) * tox * bi;
bip = bi;
bi = bim;
if (fabs (bi) > BIGNO) {
res *= BIGNI;
bi *= BIGNI;
bip *= BIGNI;
}
if (j == n) res = bip;
}
res *= (gmt_i0 (GMT, x) / bi);
if (x < 0.0 && (n%2)) res = -res;
return (res);
}
double gmt_k0 (struct GMT_CTRL *GMT, double x) {
/* Modified from code in Press et al. */
double y, z, res;
gmt_M_unused(GMT);
if (x < 0.0) x = -x;
if (x <= 2.0) {
y = 0.25 * x * x;
z = x * x / 14.0625;
res = (-log(0.5*x) * (1.0 + z * (3.5156229 + z * (3.0899424 + z * (1.2067492 + z * (0.2659732 + z * (0.360768e-1 + z * 0.45813e-2))))))) + (-M_EULER + y * (0.42278420 + y * (0.23069756 + y * (0.3488590e-1 + y * (0.262698e-2 + y * (0.10750e-3 + y * 0.74e-5))))));
}
else {
y = 2.0 / x;
res = (exp (-x) / sqrt (x)) * (1.25331414 + y * (-0.7832358e-1 + y * (0.2189568e-1 + y * (-0.1062446e-1 + y * (0.587872e-2 + y * (-0.251540e-2 + y * 0.53208e-3))))));
}
return (res);
}
double gmt_k1 (struct GMT_CTRL *GMT, double x) {
/* Modified Bessel function K1(x) */
double y, res;
if (x < 0.0) x = -x;
if (x <= 2.0) {
y = x * x / 4.0;
res = (log (0.5 * x) * gmt_i1 (GMT, x)) + (1.0 / x) * (1.0 + y * (0.15443144 + y * (-0.67278579 + y * (-0.18156897 + y * (-0.01919402 + y * (-0.00110404 - y * 0.00004686))))));
}
else {
y = 2.0 / x;
res = (exp (-x) / sqrt (x)) * (1.25331414 + y * (0.23498619 + y * (-0.03655620 + y * (0.01504268 + y * (-0.00780353 + y * (0.00325614 - y * 0.00068245))))));
}
return (res);
}
double gmt_kn (struct GMT_CTRL *GMT, unsigned int n, double x) {
/* Modified Bessel function Kn(x) */
unsigned int j;
double bkm, bk, bkp, tox;
if (n == 0) return (gmt_k0 (GMT, x));
if (n == 1) return (gmt_k1 (GMT, x));
tox = 2.0 / x;
bkm = gmt_k0 (GMT, x);
bk = gmt_k1 (GMT, x);
for (j = 1; j <= (n-1); j++) {
bkp = bkm + j * tox * bk;
bkm = bk;
bk = bkp;
}
return (bk);
}
double gmt_plm (struct GMT_CTRL *GMT, int l, int m, double x) {
/* Unnormalized associated Legendre polynomial of degree l and order m, including
* Condon-Shortley phase (-1)^m */
double fact, pll = 0, pmm, pmmp1, somx2;
int i, ll;
/* x is cosine of colatitude and must be -1 <= x <= +1 */
if (fabs(x) > 1.0) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "|x| > 1.0 in gmt_plm\n");
return (GMT->session.d_NaN);
}
if (m < 0 || m > l) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmt_plm requires 0 <= m <= l.\n");
return (GMT->session.d_NaN);
}
pmm = 1.0;
if (m > 0) {
somx2 = d_sqrt ((1.0 - x) * (1.0 + x));
fact = 1.0;
/* This loop used to go to i < m; corrected to i <= m by WHFS */
for (i = 1; i <= m; i++) {
pmm *= -fact * somx2;
fact += 2.0;
}
}
if (l == m) return (pmm);
pmmp1 = x * (2*m + 1) * pmm;
if (l == (m + 1)) return (pmmp1);
for (ll = (m+2); ll <= l; ll++) {
pll = (x * (2*ll - 1) * pmmp1 - (ll + m - 1) * pmm) / (ll - m);
pmm = pmmp1;
pmmp1 = pll;
}
return (pll);
}