Refactor verify_results and Build.PL, add Const tests for math constants like $M_PI.
[Math-GSL.git] / CDF.i
blob5dd76e9a2c375e004fb889cc09f4d761664f3885
1 %module CDF
2 %{
3 #include "/usr/local/include/gsl/gsl_cdf.h"
4 %}
6 %include "/usr/local/include/gsl/gsl_cdf.h"
8 %perlcode %{
10 our @EXPORT_OK = qw/ gsl_cdf_ugaussian_P gsl_cdf_ugaussian_Q gsl_cdf_ugaussian_Pinv
11 gsl_cdf_ugaussian_Qinv gsl_cdf_gaussian_P gsl_cdf_gaussian_Q
12 gsl_cdf_gaussian_Pinv gsl_cdf_gaussian_Qinv gsl_cdf_gamma_P
13 gsl_cdf_gamma_Q gsl_cdf_gamma_Pinv gsl_cdf_gamma_Qinv
14 gsl_cdf_cauchy_P gsl_cdf_cauchy_Q gsl_cdf_cauchy_Pinv
15 gsl_cdf_cauchy_Qinv gsl_cdf_laplace_P gsl_cdf_laplace_Q
16 gsl_cdf_laplace_Pinv gsl_cdf_laplace_Qinv gsl_cdf_rayleigh_P
17 gsl_cdf_rayleigh_Q gsl_cdf_rayleigh_Pinv gsl_cdf_rayleigh_Qinv
18 gsl_cdf_chisq_P gsl_cdf_chisq_Q gsl_cdf_chisq_Pinv
19 gsl_cdf_chisq_Qinv gsl_cdf_exponential_P gsl_cdf_exponential_Q
20 gsl_cdf_exponential_Pinv gsl_cdf_exponential_Qinv gsl_cdf_exppow_P
21 gsl_cdf_exppow_Q gsl_cdf_tdist_P gsl_cdf_tdist_Q
22 gsl_cdf_tdist_Pinv gsl_cdf_tdist_Qinv gsl_cdf_fdist_P
23 gsl_cdf_fdist_Q gsl_cdf_fdist_Pinv gsl_cdf_fdist_Qinv
24 gsl_cdf_beta_P gsl_cdf_beta_Q gsl_cdf_beta_Pinv
25 gsl_cdf_beta_Qinv gsl_cdf_flat_P gsl_cdf_flat_Q
26 gsl_cdf_flat_Pinv gsl_cdf_flat_Qinv gsl_cdf_lognormal_P
27 gsl_cdf_lognormal_Q gsl_cdf_lognormal_Pinv gsl_cdf_lognormal_Qinv
28 gsl_cdf_gumbel1_P gsl_cdf_gumbel1_Q gsl_cdf_gumbel1_Pinv
29 gsl_cdf_gumbel1_Qinv gsl_cdf_gumbel2_P gsl_cdf_gumbel2_Q
30 gsl_cdf_gumbel2_Pinv gsl_cdf_gumbel2_Qinv gsl_cdf_weibull_P
31 gsl_cdf_weibull_Q gsl_cdf_weibull_Pinv gsl_cdf_weibull_Qinv
32 gsl_cdf_pareto_P gsl_cdf_pareto_Q gsl_cdf_pareto_Pinv
33 gsl_cdf_pareto_Qinv gsl_cdf_logistic_P gsl_cdf_logistic_Q
34 gsl_cdf_logistic_Pinv gsl_cdf_logistic_Qinv gsl_cdf_binomial_P
35 gsl_cdf_binomial_Q gsl_cdf_poisson_P gsl_cdf_poisson_Q
36 gsl_cdf_geometric_P gsl_cdf_geometric_Q gsl_cdf_negative_binomial_P
37 gsl_cdf_negative_binomial_Q gsl_cdf_pascal_P gsl_cdf_pascal_Q
38 gsl_cdf_hypergeometric_P gsl_cdf_hypergeometric_Q
40 our %EXPORT_TAGS = ( all => [ @EXPORT_OK ], geometric => [ gsl_cdf_geometric_P , gsl_cdf_geometric_Q ], tdist => [ gsl_cdf_tdist_P , gsl_cdf_tdist_Q , gsl_cdf_tdist_Pinv , gsl_cdf_tdist_Qinv ], ugaussian => [ gsl_cdf_ugaussian_P , gsl_cdf_ugaussian_Q , gsl_cdf_ugaussian_Pinv , gsl_cdf_ugaussian_Qinv ], rayleigh => [ gsl_cdf_rayleigh_P , gsl_cdf_rayleigh_Q , gsl_cdf_rayleigh_Pinv , gsl_cdf_rayleigh_Qinv ], pascal => [ gsl_cdf_pascal_P , gsl_cdf_pascal_Q ], exponential => [ gsl_cdf_exponential_P , gsl_cdf_exponential_Q , gsl_cdf_exponential_Pinv , gsl_cdf_exponential_Qinv ], gumbel2 => [ gsl_cdf_gumbel2_P , gsl_cdf_gumbel2_Q , gsl_cdf_gumbel2_Pinv , gsl_cdf_gumbel2_Qinv ], gumbel1 => [ gsl_cdf_gumbel1_P , gsl_cdf_gumbel1_Q , gsl_cdf_gumbel1_Pinv , gsl_cdf_gumbel1_Qinv ], exppow => [ gsl_cdf_exppow_P , gsl_cdf_exppow_Q ], logistic => [ gsl_cdf_logistic_P , gsl_cdf_logistic_Q , gsl_cdf_logistic_Pinv , gsl_cdf_logistic_Qinv ], weibull => [ gsl_cdf_weibull_P , gsl_cdf_weibull_Q , gsl_cdf_weibull_Pinv , gsl_cdf_weibull_Qinv ], gaussian => [ gsl_cdf_gaussian_P , gsl_cdf_gaussian_Q , gsl_cdf_gaussian_Pinv , gsl_cdf_gaussian_Qinv ], poisson => [ gsl_cdf_poisson_P , gsl_cdf_poisson_Q ], beta => [ gsl_cdf_beta_P , gsl_cdf_beta_Q , gsl_cdf_beta_Pinv , gsl_cdf_beta_Qinv ], binomial => [ gsl_cdf_binomial_P , gsl_cdf_binomial_Q ], laplace => [ gsl_cdf_laplace_P , gsl_cdf_laplace_Q , gsl_cdf_laplace_Pinv , gsl_cdf_laplace_Qinv ], lognormal => [ gsl_cdf_lognormal_P , gsl_cdf_lognormal_Q , gsl_cdf_lognormal_Pinv , gsl_cdf_lognormal_Qinv ], cauchy => [ gsl_cdf_cauchy_P , gsl_cdf_cauchy_Q , gsl_cdf_cauchy_Pinv , gsl_cdf_cauchy_Qinv ], fdist => [ gsl_cdf_fdist_P , gsl_cdf_fdist_Q , gsl_cdf_fdist_Pinv , gsl_cdf_fdist_Qinv ], chisq => [ gsl_cdf_chisq_P , gsl_cdf_chisq_Q , gsl_cdf_chisq_Pinv , gsl_cdf_chisq_Qinv ], gamma => [ gsl_cdf_gamma_P , gsl_cdf_gamma_Q , gsl_cdf_gamma_Pinv , gsl_cdf_gamma_Qinv ], hypergeometric => [ gsl_cdf_hypergeometric_P , gsl_cdf_hypergeometric_Q ], negative => [ gsl_cdf_negative_binomial_P , gsl_cdf_negative_binomial_Q ], pareto => [ gsl_cdf_pareto_P , gsl_cdf_pareto_Q , gsl_cdf_pareto_Pinv , gsl_cdf_pareto_Qinv ], flat => [ gsl_cdf_flat_P , gsl_cdf_flat_Q , gsl_cdf_flat_Pinv , gsl_cdf_flat_Qinv ]);
42 __END__
44 =head1 NAME
46 Math::GSL::CDF
47 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the named distributions.
49 =head1 SYPNOPSIS
51 use Math::GSL::CDF qw /:all/;
53 =head1 DESCRIPTION
55 Here is a list of all the functions included in this module :
57 gsl_cdf_ugaussian_P($x)
58 gsl_cdf_ugaussian_Q($x)
59 gsl_cdf_ugaussian_Pinv($P)
60 gsl_cdf_ugaussian_Qinv($Q)
61 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the unit Gaussian distribution.
63 gsl_cdf_gaussian_P($x, $sigma)
64 gsl_cdf_gaussian_Q($x, $sigma)
65 gsl_cdf_gaussian_Pinv($P, $sigma)
66 gsl_cdf_gaussian_Qinv($Q, $sigma)
67 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Gaussian distribution with standard deviation $sigma.
69 gsl_cdf_gamma_P($x, $a, $b)
70 gsl_cdf_gamma_Q($x, $a, $b)
71 gsl_cdf_gamma_Pinv($P, $a, $b)
72 gsl_cdf_gamma_Qinv($Q, $a, $b)
73 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the gamma distribution with parameters $a and $b.
75 gsl_cdf_cauchy_P($x, $a)
76 gsl_cdf_cauchy_Q($x, $a)
77 gsl_cdf_cauchy_Pinv($P, $a)
78 gsl_cdf_cauchy_Qinv($Q, $a)
79 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Cauchy distribution with scale parameter $a.
81 gsl_cdf_laplace_P($x, $a)
82 gsl_cdf_laplace_Q($x, $a)
83 gsl_cdf_laplace_Pinv($P, $a)
84 gsl_cdf_laplace_Qinv($Q, $a)
85 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Laplace distribution with width $a.
87 gsl_cdf_rayleigh_P($x, $sigma)
88 gsl_cdf_rayleigh_Q($x, $sigma)
89 gsl_cdf_rayleigh_Pinv($P, $sigma)
90 gsl_cdf_rayleigh_Qinv($Q, $sigma)
91 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Rayleigh distribution with scale parameter $sigma.
93 gsl_cdf_chisq_P($x, $nu)
94 gsl_cdf_chisq_Q($x, $nu)
95 gsl_cdf_chisq_Pinv($P, $nu)
96 gsl_cdf_chisq_Qinv($Q, $nu)
97 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the chi-squared distribution with $nu degrees of freedom.
99 gsl_cdf_exponential_P($x, $mu)
100 gsl_cdf_exponential_Q($x, $mu)
101 gsl_cdf_exponential_Pinv($P, $mu)
102 gsl_cdf_exponential_Qinv($Q, $mu)
103 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Laplace distribution with width $a.
105 gsl_cdf_exppow_P($x, $a, $b)
106 gsl_cdf_exppow_Q($x, $a, $b)
107 - These functions compute the cumulative distribution functions P(x), Q(x) for the exponential power distribution with parameters $a and $b.
109 gsl_cdf_tdist_P($x, $nu)
110 gsl_cdf_tdist_Q($x, $nu)
111 gsl_cdf_tdist_Pinv($P, $nu)
112 gsl_cdf_tdist_Qinv($Q, $nu)
113 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the t-distribution with $nu degrees of freedom.
115 gsl_cdf_fdist_P($x, $nu1, $nu2)
116 gsl_cdf_fdist_Q($x, $nu1, $nu2)
117 gsl_cdf_fdist_Pinv($P, $nu1, $nu2)
118 gsl_cdf_fdist_Qinv($Q, $nu1, $nu2)
119 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the F-distribution with $nu1 and $nu2 degrees of freedom.
121 gsl_cdf_beta_P($x, $a, $b)
122 gsl_cdf_beta_Q($x, $a, $b)
123 gsl_cdf_beta_Pinv($P, $a, $b)
124 gsl_cdf_beta_Qinv($Q, $a, $b)
125 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the beta distribution with parameters $a and $b.
127 gsl_cdf_flat_P($x, $a, $b)
128 gsl_cdf_flat_Q($x, $a, $b)
129 gsl_cdf_flat_Pinv($P, $a, $b)
130 gsl_cdf_flat_Qinv($Q, $a, $b)
131 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for a uniform distribution from $a to $b.
133 gsl_cdf_lognormal_P($x, $zeta, $sigma)
134 gsl_cdf_lognormal_Q($x, $zeta, $sigma)
135 gsl_cdf_lognormal_Pinv($P, $zeta, $sigma)
136 gsl_cdf_lognormal_Qinv($Q, $zeta, $sigma)
137 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the lognormal distribution with parameters $zeta and $sigma.
139 gsl_cdf_gumbel1_P($x, $a, $b)
140 gsl_cdf_gumbel1_Q($x, $a, $b)
141 gsl_cdf_gumbel1_Pinv($P, $a, $b)
142 gsl_cdf_gumbel1_Qinv($Q, $a, $b)
143 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Type-1 Gumbel distribution with parameters $a and $b.
145 gsl_cdf_gumbel2_P($x, $a, $b)
146 gsl_cdf_gumbel2_Q($x, $a, $b)
147 gsl_cdf_gumbel2_Pinv($P, $a, $b)
148 gsl_cdf_gumbel2_Qinv($Q, $a, $b)
149 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Type-2 Gumbel distribution with parameters $a and $b.
151 gsl_cdf_weibull_P($x, $a, $b)
152 gsl_cdf_weibull_Q($x, $a, $b)
153 gsl_cdf_weibull_Pinv($P, $a, $b)
154 gsl_cdf_weibull_Qinv($Q, $a, $b)
155 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Type-1 Gumbel distribution with parameters $a and $b.
157 gsl_cdf_pareto_P($x, $a, $b)
158 gsl_cdf_pareto_Q($x, $a, $b)
159 gsl_cdf_pareto_Pinv($P, $a, $b)
160 gsl_cdf_pareto_Qinv($Q, $a, $b)
161 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the Pareto distribution with exponent $a and scale $b.
163 gsl_cdf_logistic_P($x, $a)
164 gsl_cdf_logistic_Q($x, $a)
165 gsl_cdf_logistic_Pinv($P, $a)
166 gsl_cdf_logistic_Qinv($Q, $a)
167 - These functions compute the cumulative distribution functions P(x), Q(x) and their inverses for the logistic distribution with scale parameter a.
169 gsl_cdf_binomial_P($k, $p, $n)
170 gsl_cdf_binomial_Q($k, $p, $n)
171 - These functions compute the cumulative distribution functions P(k), Q(k) for the binomial distribution with parameters $p and $n.
173 gsl_cdf_poisson_P($k, $mu)
174 gsl_cdf_poisson_Q($k, $mu)
175 - These functions compute the cumulative distribution functions P(k), Q(k) for the Poisson distribution with parameter $mu.
177 gsl_cdf_geometric_P($k, $p)
178 gsl_cdf_geometric_Q($k, $p)
179 - These functions compute the cumulative distribution functions P(k), Q(k) for the geometric distribution with parameter $p.
181 gsl_cdf_negative_binomial_P($k, $p, $n)
182 gsl_cdf_negative_binomial_Q($k, $p, $n)
183 - These functions compute the cumulative distribution functions P(k), Q(k) for the negative binomial distribution with parameters $p and $n.
185 gsl_cdf_pascal_P($k, $p, $n)
186 gsl_cdf_pascal_Q($k, $p, $n)
187 - These functions compute the cumulative distribution functions P(k), Q(k) for the Pascal distribution with parameters $p and $n.
189 gsl_cdf_hypergeometric_P($k, $n1, $n2, $t)
190 gsl_cdf_hypergeometric_Q($k, $n1, $n2, $t)
191 - These functions compute the cumulative distribution functions P(k), Q(k) for the hypergeometric distribution with parameters $n1, $n2 and $t.
194 You have to add the functions you want to use inside the qw /put_funtion_here / with spaces between each function. You can also write use Math::GSL::CDF qw/:all/ to use all avaible functions of the module. Other tags are also avaible, here is a complete list of all tags for this module :
196 geometric
197 tdist
198 ugaussian
199 rayleigh
200 pascal
201 exponential
202 gumbel2
203 gumbel1
204 exppow
205 logistic
206 weibull
207 gaussian
208 poisson
209 beta
210 binomial
211 laplace
212 lognormal
213 cauchy
214 fdist
215 chisq
216 gamma
217 hypergeometric
218 negative
219 pareto
220 flat
222 For example the beta tag contains theses functions : gsl_cdf_beta_P, gsl_cdf_beta_Q, gsl_cdf_beta_Pinv, gsl_cdf_beta_Qinv.
224 For more informations on the functions, we refer you to the GSL offcial documentation: http://www.gnu.org/software/gsl/manual/html_node/
225 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/ name_of_the_function_you_want
227 =head1 EXAMPLES
229 example using tag:
231 use Math::GSL::CDF qw /:beta/;
232 print gsl_cdf_beta_P(1,2,3) . "\n";
235 example using functions:
237 use Math::GSL::CDF qw /gsl_cdf_laplace_P gsl_cdf_laplace_Q/;
239 print gsl_cdf_laplace_P(2,3). "\n";
240 print gsl_cdf_laplace_Q(2,3). "\n";
242 =head1 AUTHOR
244 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
246 =head1 COPYRIGHT AND LICENSE
248 Copyright (C) 2008 Jonathan Leto and Thierry Moisan
250 This program is free software; you can redistribute it and/or modify it
251 under the same terms as Perl itself.
253 =cut