/
PDF.tt
189 lines (179 loc) · 6.64 KB
/
PDF.tt
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
#include "EXTERN.h"
#include "perl.h"
#include "XSUB.h"
#include "ppport.h"
#include <gsl/gsl_randist.h>
void warn_handler (const char * reason, const char * file, int line, int gsl_errno) {
warn("Error(gsl-internal): '%s/%s' (%s:%d)", gsl_strerror(gsl_errno), reason, file, line);
return;
}
size_t arr_ref_size(SV * data, char * warnmsg) {
size_t rv;
if ((!SvROK(data)) || (SvTYPE(SvRV(data)) != SVt_PVAV) || ((rv = av_len((AV *)SvRV(data))) < 0)) {
if (warnmsg) warn(warnmsg);
return -1;
}
return rv+1;
}
MODULE = Math::EasyGSL::PDF PACKAGE = Math::EasyGSL::PDF
## GSL function: double gsl_ran_dirichlet_pdf(size_t K, const double alpha[], const double theta[]);
void
pdf_dirichlet(K,alpha,theta)
size_t K;
SV * alpha;
SV * theta;
INIT:
int i;
AV * av;
double * alpha_;
size_t alpha_size;
double * theta_;
size_t theta_size;
PPCODE:
/* check valid ARRAYREF */
alpha_size = arr_ref_size(alpha, "Warning: pdf_dirichlet() - invalid 'alpha' argument");
theta_size = arr_ref_size(theta, "Warning: pdf_dirichlet() - invalid 'theta' argument");
if (alpha_size<0 || theta_size<0) XSRETURN_UNDEF;
/* check size */
if (alpha_size!=theta_size || alpha_size!=K) {
warn("Warning: pdf_dirichlet() - array sizes differ: K(%d), alpha_size(%d), theta_size(%d)", K, alpha_size, theta_size);
XSRETURN_UNDEF;
}
/* allocate memory */
alpha_ = malloc(alpha_size * sizeof(double));
theta_ = malloc(theta_size * sizeof(double));
if (!alpha_ || !theta_) {
warn("Warning: pdf_dirichlet() - malloc failed");
if (alpha_) free(alpha_);
if (theta_) free(theta_);
XSRETURN_UNDEF;
}
/* copy data */
for(i=0, av=(AV*)SvRV(alpha); i<alpha_size; i++) alpha_[i] = SvNV(*av_fetch(av,i,0));
for(i=0, av=(AV*)SvRV(theta); i<theta_size; i++) theta_[i] = SvNV(*av_fetch(av,i,0));
/* do the job */
XPUSHs(sv_2mortal(newSVnv(gsl_ran_dirichlet_pdf(K,alpha_,theta_))));
/* free memory */
free(alpha_);
free(theta_);
## GSL function: double gsl_ran_dirichlet_lnpdf(size_t K, const double alpha[], const double theta[]);
void
pdf_dirichlet_ln(K,alpha,theta)
size_t K;
SV * alpha;
SV * theta;
INIT:
int i;
AV * av;
double * alpha_;
size_t alpha_size;
double * theta_;
size_t theta_size;
PPCODE:
/* check valid ARRAYREF */
alpha_size = arr_ref_size(alpha, "Warning: pdf_dirichlet_ln() - invalid 'alpha' argument");
theta_size = arr_ref_size(theta, "Warning: pdf_dirichlet_ln() - invalid 'theta' argument");
if (alpha_size<0 || theta_size<0) XSRETURN_UNDEF;
/* check size */
if (alpha_size!=theta_size || alpha_size!=K) {
warn("Warning: pdf_dirichlet_ln() - array sizes differ: K(%d), alpha_size(%d), theta_size(%d)", K, alpha_size, theta_size);
XSRETURN_UNDEF;
}
/* allocate memory */
alpha_ = malloc(alpha_size * sizeof(double));
theta_ = malloc(theta_size * sizeof(double));
if (!alpha_ || !theta_) {
warn("Warning: pdf_dirichlet_ln() - malloc failed");
if (alpha_) free(alpha_);
if (theta_) free(theta_);
XSRETURN_UNDEF;
}
/* copy data */
for(i=0, av=(AV*)SvRV(alpha); i<alpha_size; i++) alpha_[i] = SvNV(*av_fetch(av,i,0));
for(i=0, av=(AV*)SvRV(theta); i<theta_size; i++) theta_[i] = SvNV(*av_fetch(av,i,0));
/* do the job */
XPUSHs(sv_2mortal(newSVnv(gsl_ran_dirichlet_lnpdf(K,alpha_,theta_))));
/* free memory */
free(alpha_);
free(theta_);
## GSL function: double gsl_ran_multinomial_pdf(size_t K, const double p[], const unsigned int n[]);
void
pdf_multinomial(K,p,n)
size_t K;
SV * p;
SV * n;
INIT:
int i;
AV * av;
double * p_;
size_t p_size;
unsigned int * n_;
size_t n_size;
PPCODE:
/* check valid ARRAYREF */
p_size = arr_ref_size(p, "Warning: pdf_multinomial() - invalid 'p' argument");
n_size = arr_ref_size(n, "Warning: pdf_multinomial() - invalid 'n' argument");
if (p_size<0 || n_size<0) XSRETURN_UNDEF;
/* check size */
if (p_size!=n_size || p_size!=K) {
warn("Warning: pdf_multinomial() - array sizes differ: K(%d), p_size(%d), n_size(%d)", K, p_size, n_size);
XSRETURN_UNDEF;
}
/* allocate memory */
p_ = malloc(p_size * sizeof(double));
n_ = malloc(n_size * sizeof(unsigned int));
if (!p_ || !n_) {
warn("Warning: pdf_multinomial() - malloc failed");
if (p_) free(p_);
if (n_) free(n_);
XSRETURN_UNDEF;
}
/* copy data */
for(i=0, av=(AV*)SvRV(p); i<p_size; i++) p_[i] = SvNV(*av_fetch(av,i,0));
for(i=0, av=(AV*)SvRV(n); i<n_size; i++) n_[i] = SvIV(*av_fetch(av,i,0));
/* do the job */
XPUSHs(sv_2mortal(newSVnv(gsl_ran_multinomial_pdf(K,p_,n_))));
/* free memory */
free(p_);
free(n_);
## GSL function: double gsl_ran_multinomial_lnpdf(size_t K, const double p[], const unsigned int n[]);
void
pdf_multinomial_ln(K,p,n)
size_t K;
SV * p;
SV * n;
INIT:
int i;
AV * av;
double * p_;
size_t p_size;
unsigned int * n_;
size_t n_size;
PPCODE:
/* check valid ARRAYREF */
p_size = arr_ref_size(p, "Warning: pdf_multinomial_ln() - invalid 'p' argument");
n_size = arr_ref_size(n, "Warning: pdf_multinomial_ln() - invalid 'n' argument");
if (p_size<0 || n_size<0) XSRETURN_UNDEF;
/* check size */
if (p_size!=n_size || p_size!=K) {
warn("Warning: pdf_multinomial_ln() - array sizes differ: K(%d), p_size(%d), n_size(%d)", K, p_size, n_size);
XSRETURN_UNDEF;
}
/* allocate memory */
p_ = malloc(p_size * sizeof(double));
n_ = malloc(n_size * sizeof(unsigned int));
if (!p_ || !n_) {
warn("Warning: pdf_multinomial_ln() - malloc failed");
if (p_) free(p_);
if (n_) free(n_);
XSRETURN_UNDEF;
}
/* copy data */
for(i=0, av=(AV*)SvRV(p); i<p_size; i++) p_[i] = SvNV(*av_fetch(av,i,0));
for(i=0, av=(AV*)SvRV(n); i<n_size; i++) n_[i] = SvIV(*av_fetch(av,i,0));
/* do the job */
XPUSHs(sv_2mortal(newSVnv(gsl_ran_multinomial_lnpdf(K,p_,n_))));
/* free memory */
free(p_);
free(n_);
[% INCLUDE Common.tt %]