File: Randist.pod

package info (click to toggle)
libmath-gsl-perl 0.45-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 192,156 kB
  • sloc: ansic: 895,524; perl: 24,682; makefile: 12
file content (875 lines) | stat: -rw-r--r-- 34,183 bytes parent folder | download
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
%perlcode %{

our @EXPORT_OK = qw/gsl_ran_bernoulli gsl_ran_bernoulli_pdf gsl_ran_beta
                    gsl_ran_beta_pdf gsl_ran_binomial gsl_ran_binomial_knuth
                    gsl_ran_binomial_tpe gsl_ran_binomial_pdf gsl_ran_exponential
                    gsl_ran_exponential_pdf gsl_ran_exppow gsl_ran_exppow_pdf
                    gsl_ran_cauchy gsl_ran_cauchy_pdf gsl_ran_chisq
                    gsl_ran_chisq_pdf gsl_ran_dirichlet gsl_ran_dirichlet_pdf
                    gsl_ran_dirichlet_lnpdf gsl_ran_erlang gsl_ran_erlang_pdf
                    gsl_ran_fdist gsl_ran_fdist_pdf gsl_ran_flat
                    gsl_ran_flat_pdf gsl_ran_gamma gsl_ran_gamma_int
                    gsl_ran_gamma_pdf gsl_ran_gamma_mt gsl_ran_gamma_knuth
                    gsl_ran_gaussian gsl_ran_gaussian_ratio_method gsl_ran_gaussian_ziggurat
                    gsl_ran_gaussian_pdf gsl_ran_ugaussian gsl_ran_ugaussian_ratio_method
                    gsl_ran_ugaussian_pdf gsl_ran_gaussian_tail gsl_ran_gaussian_tail_pdf
                    gsl_ran_ugaussian_tail gsl_ran_ugaussian_tail_pdf gsl_ran_bivariate_gaussian
                    gsl_ran_bivariate_gaussian_pdf gsl_ran_landau gsl_ran_landau_pdf
                    gsl_ran_geometric gsl_ran_geometric_pdf gsl_ran_hypergeometric
                    gsl_ran_hypergeometric_pdf gsl_ran_gumbel1 gsl_ran_gumbel1_pdf
                    gsl_ran_gumbel2 gsl_ran_gumbel2_pdf gsl_ran_logistic
                    gsl_ran_logistic_pdf gsl_ran_lognormal gsl_ran_lognormal_pdf
                    gsl_ran_logarithmic gsl_ran_logarithmic_pdf gsl_ran_multinomial
                    gsl_ran_multinomial_pdf gsl_ran_multinomial_lnpdf
                    gsl_ran_negative_binomial gsl_ran_negative_binomial_pdf gsl_ran_pascal
                    gsl_ran_pascal_pdf gsl_ran_pareto gsl_ran_pareto_pdf
                    gsl_ran_poisson gsl_ran_poisson_array
                    gsl_ran_poisson_pdf gsl_ran_rayleigh gsl_ran_rayleigh_pdf
                    gsl_ran_rayleigh_tail gsl_ran_rayleigh_tail_pdf gsl_ran_tdist
                    gsl_ran_tdist_pdf gsl_ran_laplace gsl_ran_laplace_pdf
                    gsl_ran_levy gsl_ran_levy_skew gsl_ran_weibull
                    gsl_ran_weibull_pdf gsl_ran_dir_2d gsl_ran_dir_2d_trig_method
                    gsl_ran_dir_3d gsl_ran_dir_nd
                    gsl_ran_discrete_t gsl_ran_discrete_free gsl_ran_discrete_preproc
                    gsl_ran_discrete gsl_ran_discrete_pdf
                /;

our %EXPORT_TAGS = (
        all                 => [  @EXPORT_OK ],
        logarithmic         => [ gsl_ran_logarithmic , gsl_ran_logarithmic_pdf ],
        exponential         => [ gsl_ran_exponential , gsl_ran_exponential_pdf ],
        gumbel1             => [ gsl_ran_gumbel1 , gsl_ran_gumbel1_pdf ],
        exppow              => [ gsl_ran_exppow , gsl_ran_exppow_pdf ],
        logistic            => [ gsl_ran_logistic , gsl_ran_logistic_pdf ],
        gaussian            => [
                                 gsl_ran_gaussian , gsl_ran_gaussian_ratio_method , gsl_ran_gaussian_ziggurat,
                                 gsl_ran_gaussian_pdf , gsl_ran_gaussian_tail , gsl_ran_gaussian_tail_pdf
                               ],
        poisson             => [ gsl_ran_poisson , gsl_ran_poisson_array , gsl_ran_poisson_pdf ],
        binomial            => [ gsl_ran_binomial , gsl_ran_binomial_knuth , gsl_ran_binomial_tpe ,
                                 gsl_ran_binomial_pdf ],
        fdist               => [ gsl_ran_fdist , gsl_ran_fdist_pdf ],
        chisq               => [ gsl_ran_chisq , gsl_ran_chisq_pdf ],
        gamma               => [ gsl_ran_gamma , gsl_ran_gamma_int , gsl_ran_gamma_pdf , gsl_ran_gamma_mt , gsl_ran_gamma_knuth ],
        hypergeometric      => [ gsl_ran_hypergeometric , gsl_ran_hypergeometric_pdf ],
        dirichlet           => [ gsl_ran_dirichlet , gsl_ran_dirichlet_pdf , gsl_ran_dirichlet_lnpdf ],
        negative            => [ gsl_ran_negative_binomial , gsl_ran_negative_binomial_pdf ],
        flat                => [ gsl_ran_flat , gsl_ran_flat_pdf ],
        geometric           => [ gsl_ran_geometric , gsl_ran_geometric_pdf ],
        discrete            => [ gsl_ran_discrete , gsl_ran_discrete_pdf],
        tdist               => [ gsl_ran_tdist , gsl_ran_tdist_pdf ],
        ugaussian           => [ gsl_ran_ugaussian , gsl_ran_ugaussian_ratio_method , gsl_ran_ugaussian_pdf ,
                                 gsl_ran_ugaussian_tail , gsl_ran_ugaussian_tail_pdf ],
        rayleigh            => [ gsl_ran_rayleigh , gsl_ran_rayleigh_pdf , gsl_ran_rayleigh_tail ,
                                 gsl_ran_rayleigh_tail_pdf ],
        dir                 => [ gsl_ran_dir_2d , gsl_ran_dir_2d_trig_method , gsl_ran_dir_3d , gsl_ran_dir_nd ],
        pascal              => [ gsl_ran_pascal , gsl_ran_pascal_pdf ],
        gumbel2             => [ gsl_ran_gumbel2 , gsl_ran_gumbel2_pdf ],
        landau              => [ gsl_ran_landau , gsl_ran_landau_pdf ],
        bernoulli           => [ gsl_ran_bernoulli , gsl_ran_bernoulli_pdf ],
        weibull             => [ gsl_ran_weibull , gsl_ran_weibull_pdf ],
        multinomial         => [ gsl_ran_multinomial , gsl_ran_multinomial_pdf , gsl_ran_multinomial_lnpdf ],
        beta                => [ gsl_ran_beta , gsl_ran_beta_pdf ],
        lognormal           => [ gsl_ran_lognormal , gsl_ran_lognormal_pdf ],
        laplace             => [ gsl_ran_laplace , gsl_ran_laplace_pdf ],
        erlang              => [ gsl_ran_erlang , gsl_ran_erlang_pdf ],
        cauchy              => [ gsl_ran_cauchy , gsl_ran_cauchy_pdf ],
        levy                => [ gsl_ran_levy , gsl_ran_levy_skew ],
        bivariate           => [ gsl_ran_bivariate_gaussian , gsl_ran_bivariate_gaussian_pdf ],
        pareto              => [ gsl_ran_pareto , gsl_ran_pareto_pdf ]
);

__END__

=encoding utf8

=head1 NAME

Math::GSL::Randist - Probability Distributions

=head1 SYNOPSIS

 use Math::GSL::RNG;
 use Math::GSL::Randist qw/:all/;

 my $rng = Math::GSL::RNG->new();
 my $coinflip = gsl_ran_bernoulli($rng->raw(), .5);


=head1 DESCRIPTION

Here is a list of all the functions included in this module. For all sampling methods, the first argument $r is a raw gsl_rnd structure retrieve by calling raw() on an Math::GSL::RNG object.

=head2 Bernoulli

=over

=item gsl_ran_bernoulli($r, $p)

This function returns either 0 or 1, the result of a Bernoulli trial with probability $p. The probability distribution for a Bernoulli trial is, p(0) = 1 - $p and  p(1) = $p. $r is a gsl_rng structure.

=item gsl_ran_bernoulli_pdf($k, $p)

This function computes the probability p($k) of obtaining $k from a Bernoulli distribution with probability parameter $p, using the formula given above.

=back

=head2 Beta

=over

=item gsl_ran_beta($r, $a, $b)

This function returns a random variate from the beta distribution. The distribution function is, p($x) dx = {Gamma($a+$b) \ Gamma($a) Gamma($b)} $x**{$a-1} (1-$x)**{$b-1} dx for 0 <= $x <= 1.$r is a gsl_rng structure.

=item gsl_ran_beta_pdf($x, $a, $b)

This function computes the probability density p($x) at $x for a beta distribution with parameters $a and $b, using the formula given above.

=back

=head2 Binomial

=over

=item gsl_ran_binomial($k, $p, $n)

This function returns a random integer from the binomial distribution, the number of successes in n independent trials with probability $p. The probability distribution for binomial variates is p($k) = {$n! \ $k! ($n-$k)! } $p**$k (1-$p)^{$n-$k} for 0 <= $k <= $n.
Uses Binomial Triangle Parallelogram Exponential algorithm.

=item gsl_ran_binomial_knuth($k, $p, $n)

Alternative and original implementation for gsl_ran_binomial using Knuth's algorithm.

=item gsl_ran_binomial_tpe($k, $p, $n)

Same as gsl_ran_binomial.

=item gsl_ran_binomial_pdf($k, $p, $n)

This function computes the probability p($k) of obtaining $k from a binomial distribution with parameters $p and $n, using the formula given above.

=back

=head2 Exponential

=over

=item gsl_ran_exponential($r, $mu)

This function returns a random variate from the exponential distribution with mean $mu. The distribution is, p($x) dx = {1 \ $mu} exp(-$x/$mu) dx for $x >= 0. $r is a gsl_rng structure.

=item gsl_ran_exponential_pdf($x, $mu)

This function computes the probability density p($x) at $x for an exponential distribution with mean $mu, using the formula given above.

=back

=head2 Exponential Power

=over

=item gsl_ran_exppow($r, $a, $b)

This function returns a random variate from the exponential power distribution with scale parameter $a and exponent $b. The distribution is, p(x) dx = {1 / 2 $a Gamma(1+1/$b)} exp(-|$x/$a|**$b) dx for $x >= 0. For $b = 1 this reduces to the Laplace distribution. For $b = 2 it has the same form as a gaussian distribution, but with $a = sqrt(2) sigma. $r is a gsl_rng structure.

=item gsl_ran_exppow_pdf($x, $a, $b)

This function computes the probability density p($x) at $x for an exponential power distribution with scale parameter $a and exponent $b, using the formula given above.

=back

=head2 Cauchy

=over

=item gsl_ran_cauchy($r, $scale)

This function returns a random variate from the Cauchy distribution with
$scale. The probability distribution for Cauchy random variates is,

 p(x) dx = {1 / $scale pi (1 + (x/$$scale)**2) } dx

for x in the range -infinity to +infinity.  The Cauchy distribution is also
known as the Lorentz distribution. $r is a gsl_rng structure.

=item gsl_ran_cauchy_pdf($x, $scale)

This function computes the probability density p($x) at $x for a Cauchy
distribution with $scale, using the formula given above.

=back

=head2 Chi-Squared

=over

=item gsl_ran_chisq($r, $nu)

This function returns a random variate from the chi-squared distribution with $nu degrees of freedom. The distribution function is, p(x) dx = {1 / 2 Gamma($nu/2) } (x/2)**{$nu/2 - 1} exp(-x/2) dx for $x >= 0. $r is a gsl_rng structure.

=item gsl_ran_chisq_pdf($x, $nu)

This function computes the probability density p($x) at $x for a chi-squared distribution with $nu degrees of freedom, using the formula given above.

=back

=head2 Dirichlet

=over

=item gsl_ran_dirichlet($r, $alpha)

This function returns an array of K (where K = length of $alpha array) random
variates from a Dirichlet distribution of order K-1. The distribution function
is

  p(\theta_1, ..., \theta_K) d\theta_1 ... d\theta_K =
     (1/Z) \prod_{i=1}^K \theta_i^{\alpha_i - 1} \delta(1 -\sum_{i=1}^K \theta_i) d\theta_1 ... d\theta_K

for theta_i >= 0 and alpha_i > 0. The delta function ensures that \sum \theta_i
= 1. The normalization factor Z is

  Z = {\prod_{i=1}^K \Gamma(\alpha_i)} / {\Gamma( \sum_{i=1}^K \alpha_i)}

The random variates are generated by sampling K values from gamma distributions
with parameters a=alpha_i, b=1, and renormalizing. See A.M. Law, W.D. Kelton,
Simulation Modeling and Analysis (1991).

=item gsl_ran_dirichlet_pdf($theta, $alpha)

This function computes the probability density p(\theta_1, ... , \theta_K) at
theta[K] for a Dirichlet distribution with parameters alpha[K], using the
formula given above. $alpha and $theta should be array references of the same size.
Theta should be normalized to sum to 1.

=item gsl_ran_dirichlet_lnpdf($theta, $alpha)

This function computes the logarithm of the probability density p(\theta_1, ...
, \theta_K) for a Dirichlet distribution with parameters alpha[K]. $alpha and
$theta should be array references of the same size.
Theta should be normalized to sum to 1.

=back

=head2 Erlang

=over

=item gsl_ran_erlang($r, $scale, $shape)

Equivalent to gsl_ran_gamma($r, $shape, $scale) where $shape is an integer.

=item gsl_ran_erlang_pdf

Equivalent to gsl_ran_gamma_pdf($r, $shape, $scale) where $shape is an integer.

=for comment
This is undoc'd in the main GSL manual for some reason...

=back

=head2 F-distribution

=over

=item gsl_ran_fdist($r, $nu1, $nu2)

This function returns a random variate from the F-distribution with degrees of freedom nu1 and nu2. The distribution function is, p(x) dx = { Gamma(($nu_1 + $nu_2)/2) / Gamma($nu_1/2) Gamma($nu_2/2) } $nu_1**{$nu_1/2} $nu_2**{$nu_2/2} x**{$nu_1/2 - 1} ($nu_2 + $nu_1 x)**{-$nu_1/2 -$nu_2/2} for $x >= 0. $r is a gsl_rng structure.

=item gsl_ran_fdist_pdf($x, $nu1, $nu2)

This function computes the probability density p(x) at x for an F-distribution with nu1 and nu2 degrees of freedom, using the formula given above.

=back

=head2 Uniform/Flat distribution

=over

=item gsl_ran_flat($r, $a, $b)

This function returns a random variate from the flat (uniform) distribution from a to b. The distribution is, p(x) dx = {1 / ($b-$a)} dx if $a <= x < $b and 0 otherwise. $r is a gsl_rng structure.

=item gsl_ran_flat_pdf($x, $a, $b)

This function computes the probability density p($x) at $x for a uniform distribution from $a to $b, using the formula given above.

=back

=head2 Gamma

=over

=item gsl_ran_gamma($r, $shape, $scale)

This function returns a random variate from the gamma distribution. The distribution function is,
          p(x) dx = {1 \over \Gamma($shape) $scale^$shape} x^{$shape-1} e^{-x/$scale} dx
for x > 0.
Uses Marsaglia-Tsang method. Can also be called as gsl_ran_gamma_mt.

=item gsl_ran_gamma_pdf($x, $shape, $scale)

This function computes the probability density p($x) at $x for a gamma distribution with parameters $shape and $scale, using the formula given above.

=item gsl_ran_gamma($r, $shape, $scale)

Same as gsl_ran_gamma.

=item gsl_ran_gamma_knuth($r, $shape, $scale)

Alternative implementation for gsl_ran_gamma, using algorithm in Knuth volume 2.

=for comment
gsl_ran_gamma_int, gsl_ran_gamma_mt removed because those were just internal functions for gsl_ran_gamma

=back

=head2 Gaussian/Normal

=over

=item gsl_ran_gaussian($r, $sigma)

This function returns a Gaussian random variate, with mean zero and standard deviation $sigma. The probability distribution for Gaussian random variates is, p(x) dx = {1 / sqrt{2 pi $sigma**2}} exp(-x**2 / 2 $sigma**2) dx for x in the range -infinity to +infinity. $r is a gsl_rng structure.
Uses Box-Mueller (polar) method.

=item gsl_ran_gaussian_ratio_method($r, $sigma)

This function computes a Gaussian random variate using the alternative Kinderman-Monahan-Leva ratio method.

=item gsl_ran_gaussian_ziggurat($r, $sigma)

This function computes a Gaussian random variate using the alternative Marsaglia-Tsang ziggurat ratio method. The Ziggurat algorithm is the fastest available algorithm in most cases. $r is a gsl_rng structure.

=item gsl_ran_gaussian_pdf($x, $sigma)

This function computes the probability density p($x) at $x for a Gaussian distribution with standard deviation sigma, using the formula given above.

=item gsl_ran_ugaussian($r)

=item gsl_ran_ugaussian_ratio_method($r)

=item gsl_ran_ugaussian_pdf($x)

This function computes results for the unit Gaussian distribution. It is equivalent to the gaussian functions above with a standard deviation of one, sigma = 1.

=item gsl_ran_bivariate_gaussian($r, $sigma_x, $sigma_y, $rho)

This function generates a pair of correlated Gaussian variates, with mean zero, correlation coefficient rho and standard deviations $sigma_x and $sigma_y in the x and y directions. The first value returned is x and the second y. The probability distribution for bivariate Gaussian random variates is, p(x,y) dx dy = {1 / 2 pi $sigma_x $sigma_y sqrt{1-$rho**2}} exp (-(x**2/$sigma_x**2 + y**2/$sigma_y**2 - 2 $rho x y/($sigma_x $sigma_y))/2(1- $rho**2)) dx dy for x,y in the range -infinity to +infinity. The correlation coefficient $rho should lie between 1 and -1. $r is a gsl_rng structure.

=item gsl_ran_bivariate_gaussian_pdf($x, $y, $sigma_x, $sigma_y, $rho)

This function computes the probability density p($x,$y) at ($x,$y) for a bivariate Gaussian distribution with standard deviations $sigma_x, $sigma_y and correlation coefficient $rho, using the formula given above.

=back

=head2 Gaussian Tail

=over

=item gsl_ran_gaussian_tail($r, $a, $sigma)

This function provides random variates from the upper tail of a Gaussian distribution with standard deviation sigma. The values returned are larger than the lower limit a, which must be positive. The probability distribution for Gaussian tail random variates is, p(x) dx = {1 / N($a; $sigma) sqrt{2 pi sigma**2}} exp(- x**2/(2 sigma**2)) dx for x > $a where N($a; $sigma) is the normalization constant, N($a; $sigma) = (1/2) erfc($a / sqrt(2 $sigma**2)). $r is a gsl_rng structure.

=item gsl_ran_gaussian_tail_pdf($x, $a, $gaussian)

This function computes the probability density p($x) at $x for a Gaussian tail distribution with standard deviation sigma and lower limit $a, using the formula given above.

=item gsl_ran_ugaussian_tail($r, $a)

This functions compute results for the tail of a unit Gaussian distribution. It is equivalent to the functions above with a standard deviation of one, $sigma = 1. $r is a gsl_rng structure.

=item gsl_ran_ugaussian_tail_pdf($x, $a)

This functions compute results for the tail of a unit Gaussian distribution. It is equivalent to the functions above with a standard deviation of one, $sigma = 1.

=back

=head2 Landau

=over

=item gsl_ran_landau($r)

This function returns a random variate from the Landau distribution. The probability distribution for Landau random variates is defined analytically by the complex integral, p(x) = (1/(2 \pi i)) \int_{c-i\infty}^{c+i\infty} ds exp(s log(s) + x s) For numerical purposes it is more convenient to use the following equivalent form of the integral, p(x) = (1/\pi) \int_0^\infty dt \exp(-t \log(t) - x t) \sin(\pi t). $r is a gsl_rng structure.

=item gsl_ran_landau_pdf($x)

This function computes the probability density p($x) at $x for the Landau distribution using an approximation to the formula given above.

=back

=head2 Geometric

=over

=item gsl_ran_geometric($r, $p)

This function returns a random integer from the geometric distribution, the number of independent trials with probability $p until the first success. The probability distribution for geometric variates is, p(k) =  p (1-$p)^(k-1) for k >= 1. Note that the distribution begins with k=1 with this definition. There is another convention in which the exponent k-1 is replaced by k. $r is a gsl_rng structure.

=item gsl_ran_geometric_pdf($k, $p)

This function computes the probability p($k) of obtaining $k from a geometric distribution with probability parameter p, using the formula given above.

=back

=head2 Hypergeometric

=over

=item gsl_ran_hypergeometric($r, $n1, $n2, $t)

This function returns a random integer from the hypergeometric distribution. The probability distribution for hypergeometric random variates is, p(k) =  C(n_1, k) C(n_2, t - k) / C(n_1 + n_2, t) where C(a,b) = a!/(b!(a-b)!) and t <= n_1 + n_2. The domain of k is max(0,t-n_2), ..., min(t,n_1). If a population contains n_1 elements of "type 1" and n_2 elements of "type 2" then the hypergeometric distribution gives the probability of obtaining k elements of "type 1" in t samples from the population without replacement. $r is a gsl_rng structure.

=item gsl_ran_hypergeometric_pdf($k, $n1, $n2, $t)

This function computes the probability p(k) of obtaining k from a hypergeometric distribution with parameters $n1, $n2 $t, using the formula given above.

=back

=head2 Gumbel

=over

=item gsl_ran_gumbel1($r, $a, $b)

This function returns a random variate from the Type-1 Gumbel distribution. The Type-1 Gumbel distribution function is, p(x) dx = a b \exp(-(b \exp(-ax) + ax)) dx for -\infty < x < \infty. $r is a gsl_rng structure.

=item gsl_ran_gumbel1_pdf($x, $a, $b)

This function computes the probability density p($x) at $x for a Type-1 Gumbel distribution with parameters $a and $b, using the formula given above.

=item gsl_ran_gumbel2($r, $a, $b)

This function returns a random variate from the Type-2 Gumbel distribution. The Type-2 Gumbel distribution function is, p(x) dx = a b x^{-a-1} \exp(-b x^{-a}) dx for 0 < x < \infty. $r is a gsl_rng structure.

=item gsl_ran_gumbel2_pdf($x, $a, $b)

This function computes the probability density p($x) at $x for a Type-2 Gumbel distribution with parameters $a and $b, using the formula given above.

=back

=head2 Logistic

=over

=item gsl_ran_logistic($r, $a)

This function returns a random variate from the logistic distribution. The distribution function is, p(x) dx = { \exp(-x/a) \over a (1 + \exp(-x/a))^2 } dx for -\infty < x < +\infty. $r is a gsl_rng structure.

=item gsl_ran_logistic_pdf($x, $a)

This function computes the probability density p($x) at $x for a logistic distribution with scale parameter $a, using the formula given above.

=back

=head2 Lognormal

=over

=item gsl_ran_lognormal($r, $zeta, $sigma)

This function returns a random variate from the lognormal distribution. The distribution function is, p(x) dx = {1 \over x \sqrt{2 \pi \sigma^2} } \exp(-(\ln(x) - \zeta)^2/2 \sigma^2) dx for x > 0. $r is a gsl_rng structure.

=item gsl_ran_lognormal_pdf($x, $zeta, $sigma)

This function computes the probability density p($x) at $x for a lognormal distribution with parameters $zeta and $sigma, using the formula given above.

=back

=head2 Logarithmic

=over

=item gsl_ran_logarithmic($r, $p)

This function returns a random integer from the logarithmic distribution. The probability distribution for logarithmic random variates is, p(k) = {-1 \over \log(1-p)} {(p^k \over k)} for k >= 1. $r is a gsl_rng structure.

=item gsl_ran_logarithmic_pdf($k, $p)

This function computes the probability p($k) of obtaining $k from a logarithmic distribution with probability parameter $p, using the formula given above.

=back

=head2 Multinomial

=over

=item gsl_ran_multinomial($r, $P, $N)

This function computes and returns a random sample n[] from the multinomial
distribution formed by N trials from an underlying distribution p[K]. The
distribution function for n[] is,

 P(n_1, n_2, ..., n_K) =
    (N!/(n_1! n_2! ... n_K!)) p_1^n_1 p_2^n_2 ... p_K^n_K

where (n_1, n_2, ..., n_K) are nonnegative integers with sum_{k=1}^K n_k = N,
and (p_1, p_2, ..., p_K) is a probability distribution with \sum p_i = 1. If
the array p[K] is not normalized then its entries will be treated as weights
and normalized appropriately.

Random variates are generated using the conditional binomial method (see C.S.
Davis, The computer generation of multinomial random variates, Comp. Stat. Data
Anal. 16 (1993) 205-217 for details).

=item gsl_ran_multinomial_pdf($counts, $P)

This function returns the probability for the multinomial
distribution P(counts[1], counts[2], ..., counts[K]) with parameters p[K].

=item gsl_ran_multinomial_lnpdf($counts, $P)

This function returns the logarithm of the probability for the multinomial
distribution P(counts[1], counts[2], ..., counts[K]) with parameters p[K].

=back

=head2 Negative Binomial

=over

=item gsl_ran_negative_binomial($r, $p, $n)

This function returns a random integer from the negative binomial distribution, the number of failures occurring before n successes in independent trials with probability p of success. The probability distribution for negative binomial variates is, p(k) = {\Gamma(n + k) \over \Gamma(k+1) \Gamma(n) } p^n (1-p)^k Note that n is not required to be an integer.

=item gsl_ran_negative_binomial_pdf($k, $p, $n)

This function computes the probability p($k) of obtaining $k from a negative binomial distribution with parameters $p and $n, using the formula given above.

=back

=head2 Pascal

=over

=item gsl_ran_pascal($r, $p, $n)

This function returns a random integer from the Pascal distribution. The Pascal distribution is simply a negative binomial distribution with an integer value of $n. p($k) = {($n + $k - 1)! \ $k! ($n - 1)! } $p**$n (1-$p)**$k for $k >= 0. $r is gsl_rng structure

=item gsl_ran_pascal_pdf($k, $p, $n)

This function computes the probability p($k) of obtaining $k from a Pascal distribution with parameters $p and $n, using the formula given above.

=back

=head2 Pareto

=over

=item gsl_ran_pareto($r, $a, $b)

This function returns a random variate from the Pareto distribution of order $a. The distribution function is p($x) dx = ($a/$b) / ($x/$b)^{$a+1} dx for $x >= $b. $r is a gsl_rng structure

=item gsl_ran_pareto_pdf($x, $a, $b)

This function computes the probability density p(x) at x for a Pareto distribution with exponent a and scale b, using the formula given above.

=back

=head2 Poisson

=over

=item gsl_ran_poisson($r, $lambda)

This function returns a random integer from the Poisson distribution with mean
$lambda. $r is a gsl_rng structure. The probability distribution for Poisson
variates is,

 p(k) = {$lambda**$k \ $k!} exp(-$lambda)

for $k >= 0. $r is a gsl_rng structure.

=item gsl_ran_poisson_pdf($k, $lambda)

This function computes the probability p($k) of obtaining $k from a Poisson
distribution with mean $lambda, using the formula given above.

=for Comment
gsl_ran_poisson_array removed b/c it's an undoc'd function which merely returns array of poisson draws

=back

=head2 Rayleigh

=over

=item gsl_ran_rayleigh($r, $sigma)

This function returns a random variate from the Rayleigh distribution with scale parameter sigma. The distribution is, p(x) dx = {x \over \sigma^2} \exp(- x^2/(2 \sigma^2)) dx for x > 0. $r is a gsl_rng structure

=item gsl_ran_rayleigh_pdf($x, $sigma)

This function computes the probability density p($x) at $x for a Rayleigh distribution with scale parameter sigma, using the formula given above.

=item gsl_ran_rayleigh_tail($r, $a, $sigma)

This function returns a random variate from the tail of the Rayleigh distribution with scale parameter $sigma and a lower limit of $a. The distribution is, p(x) dx = {x \over \sigma^2} \exp ((a^2 - x^2) /(2 \sigma^2)) dx for x > a. $r is a gsl_rng structure

=item gsl_ran_rayleigh_tail_pdf($x, $a, $sigma)

This function computes the probability density p($x) at $x for a Rayleigh tail distribution with scale parameter $sigma and lower limit $a, using the formula given above.

=back

=head2 Student-t

=over

=item gsl_ran_tdist($r, $nu)

This function returns a random variate from the t-distribution. The distribution function is, p(x) dx = {\Gamma((\nu + 1)/2) \over \sqrt{\pi \nu} \Gamma(\nu/2)} (1 + x^2/\nu)^{-(\nu + 1)/2} dx for -\infty < x < +\infty.

=item gsl_ran_tdist_pdf($x, $nu)

This function computes the probability density p($x) at $x for a t-distribution with nu degrees of freedom, using the formula given above.

=back

=head2 Laplace

=over

=item gsl_ran_laplace($r, $a)

This function returns a random variate from the Laplace distribution with width $a. The distribution is, p(x) dx = {1 \over 2 a}  \exp(-|x/a|) dx for -\infty < x < \infty.

=item gsl_ran_laplace_pdf($x, $a)

This function computes the probability density p($x) at $x for a Laplace distribution with width $a, using the formula given above.

=back

=head2 Levy

=over

=item gsl_ran_levy($r, $c, $alpha)

This function returns a random variate from the Levy symmetric stable distribution with scale $c and exponent $alpha. The symmetric stable probability distribution is defined by a fourier transform, p(x) = {1 \over 2 \pi} \int_{-\infty}^{+\infty} dt \exp(-it x - |c t|^alpha) There is no explicit solution for the form of p(x) and the library does not define a corresponding pdf function. For \alpha = 1 the distribution reduces to the Cauchy distribution. For \alpha = 2 it is a Gaussian distribution with \sigma = \sqrt{2} c. For \alpha < 1 the tails of the distribution become extremely wide. The algorithm only works for 0 < alpha <= 2. $r is a gsl_rng structure

=item gsl_ran_levy_skew($r, $c, $alpha, $beta)

This function returns a random variate from the Levy skew stable distribution with scale $c, exponent $alpha and skewness parameter $beta. The skewness parameter must lie in the range [-1,1]. The Levy skew stable probability distribution is defined by a fourier transform, p(x) = {1 \over 2 \pi} \int_{-\infty}^{+\infty} dt \exp(-it x - |c t|^alpha (1-i beta sign(t) tan(pi alpha/2))) When \alpha = 1 the term \tan(\pi \alpha/2) is replaced by -(2/\pi)\log|t|. There is no explicit solution for the form of p(x) and the library does not define a corresponding pdf function. For $alpha = 2 the distribution reduces to a Gaussian distribution with $sigma = sqrt(2) $c and the skewness parameter has no effect. For $alpha < 1 the tails of the distribution become extremely wide. The symmetric distribution corresponds to $beta = 0. The algorithm only works for 0 < $alpha <= 2. The Levy alpha-stable distributions have the property that if N alpha-stable variates are drawn from the distribution p(c, \alpha, \beta) then the sum Y = X_1 + X_2 + \dots + X_N will also be distributed as an alpha-stable variate, p(N^(1/\alpha) c, \alpha, \beta). $r is a gsl_rng structure

=back

=head2 Weibull

=over

=item gsl_ran_weibull($r, $scale, $exponent)

This function returns a random variate from the Weibull distribution with
$scale and $exponent (aka scale). The distribution function is

 p(x) dx = {$exponent \over $scale^$exponent} x^{$exponent-1}
           \exp(-(x/$scale)^$exponent) dx

for x >= 0. $r is a gsl_rng structure

=item gsl_ran_weibull_pdf($x, $scale, $exponent)

This function computes the probability density p($x) at $x for a Weibull
distribution with $scale and $exponent, using the formula given above.

=back

=head2 Spherical Vector

=over

=item gsl_ran_dir_2d($r)

This function returns two values. The first is $x and the second is $y of a
random direction vector v = ($x,$y) in two dimensions. The vector is normalized
such that |v|^2 = $x^2 + $y^2 = 1. $r is a gsl_rng structure

=item gsl_ran_dir_2d_trig_method($r)

This function returns two values. The first is $x and the second is $y of a
random direction vector v = ($x,$y) in two dimensions. The vector is normalized
such that |v|^2 = $x^2 + $y^2 = 1. $r is a gsl_rng structure

=item gsl_ran_dir_3d($r)

This function returns three values. The first is $x, the second $y and the third
$z of a random direction vector v = ($x,$y,$z) in three dimensions. The vector
is normalized such that |v|^2 = x^2 + y^2 + z^2 = 1. The method employed is due
to Robert E. Knop (CACM 13, 326 (1970)), and explained in Knuth, v2, 3rd ed,
p136. It uses the surprising fact that the distribution projected along any axis
is actually uniform (this is only true for 3 dimensions).

=item gsl_ran_dir_nd (Not yet implemented )

This function returns a random direction vector v = (x_1,x_2,...,x_n) in n
dimensions. The vector is normalized such that

    |v|^2 = x_1^2 + x_2^2 + ... + x_n^2 = 1.

The method uses the fact that a multivariate Gaussian distribution is
spherically symmetric. Each component is generated to have a Gaussian
distribution, and then the components are normalized. The method is described by
Knuth, v2, 3rd ed, p135-136, and attributed to G. W. Brown, Modern Mathematics
for the Engineer (1956).

=back

=head2 Shuffling and Sampling

=over

=item gsl_ran_shuffle

Please use the C<shuffle> method in the L<GSL::RNG> module OO interface.

=item gsl_ran_choose

Please use the C<choose> method in the L<GSL::RNG> module OO interface.

=item gsl_ran_sample

Please use the C<sample> method in the L<GSL::RNG> module OO interface.

=item gsl_ran_discrete_preproc

=item gsl_ran_discrete($r, $g)

After gsl_ran_discrete_preproc has been called, you use this function to get the discrete random numbers. $r is a gsl_rng structure and $g is a gsl_ran_discrete structure

=item gsl_ran_discrete_pdf($k, $g)

Returns the probability P[$k] of observing the variable $k. Since P[$k] is not stored as part of the lookup table, it must be recomputed; this computation takes O(K), so if K is large and you care about the original array P[$k] used to create the lookup table, then you should just keep this original array P[$k] around. $r is a gsl_rng structure and $g is a gsl_ran_discrete structure

=item gsl_ran_discrete_free($g)

De-allocates the gsl_ran_discrete pointed to by g.

=back

 You have to add the functions you want to use inside the qw /put_function_here /.
 You can also write use Math::GSL::Randist qw/:all/; to use all available functions of the module.
 Other tags are also available, here is a complete list of all tags for this module :

=over

=item logarithmic

=item choose

=item exponential

=item gumbel1

=item exppow

=item sample

=item logistic

=item gaussian

=item poisson

=item binomial

=item fdist

=item chisq

=item gamma

=item hypergeometric

=item dirichlet

=item negative

=item flat

=item geometric

=item discrete

=item tdist

=item ugaussian

=item rayleigh

=item dir

=item pascal

=item gumbel2

=item shuffle

=item landau

=item bernoulli

=item weibull

=item multinomial

=item beta

=item lognormal

=item laplace

=item erlang

=item cauchy

=item levy

=item bivariate

=item pareto

=back

 For example the beta tag contains theses functions : gsl_ran_beta, gsl_ran_beta_pdf.

For more information on the functions, we refer you to the GSL official documentation:
L<http://www.gnu.org/software/gsl/manual/html_node/>


 You might also want to write

    use Math::GSL::RNG qw/:all/;

since a lot of the functions of Math::GSL::Randist take as argument a structure
that is created by Math::GSL::RNG.  Refer to Math::GSL::RNG documentation to
see how to create such a structure.

Math::GSL::CDF also contains a structure named gsl_ran_discrete_t. An example is
given in the EXAMPLES part on how to use the function related to this structure.


=head1 EXAMPLES

    use Math::GSL::Randist qw/:all/;
    print gsl_ran_exponential_pdf(5,2) . "\n";

    use Math::GSL::Randist qw/:all/;
    my $x = Math::GSL::gsl_ran_discrete_t::new;

=head1 AUTHORS

Jonathan "Duke" Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>

=head1 COPYRIGHT AND LICENSE

Copyright (C) 2008-2024 Jonathan "Duke" Leto and Thierry Moisan

This program is free software; you can redistribute it and/or modify it
under the same terms as Perl itself.

=cut
%}