File Coverage

blib/lib/Statistics/Multtest.pm
Criterion Covered Total %
statement 187 197 94.9
branch 26 38 68.4
condition 12 24 50.0
subroutine 32 32 100.0
pod 0 6 0.0
total 257 297 86.5


line stmt bran cond sub pod time code
1             package Statistics::Multtest;
2              
3 6     6   233318 use List::Vectorize;
  6         73802  
  6         1772  
4 6     6   55 use Carp;
  6         17  
  6         386  
5 6     6   38 use strict;
  6         12  
  6         181  
6 6     6   58 use constant {EPS => 1e-10};
  6         23  
  6         941  
7              
8             require Exporter;
9             our @ISA = qw(Exporter);
10             our @EXPORT_OK = qw(bonferroni holm hommel hochberg BH BY qvalue);
11             our %EXPORT_TAGS = (all => [qw(bonferroni holm hommel hochberg BH BY qvalue)]);
12              
13             our $VERSION = '0.14';
14              
15             1;
16              
17             BEGIN {
18 6     6   44 no strict 'refs';
  6         13  
  6         836  
19 6     6   29 for my $adjp (qw(bonferroni holm hommel hochberg BH BY qvalue)) {
20 42         12771 *{$adjp} = sub {
21 20     20   40037 my $p = shift;
22 20         45 my $name;
23 20         69 ($name, $p) = initial($p);
24 20         38 *private = *{"_$adjp"};
  20         157  
25 20         89 my $adjp = private($p, @_);
26 18         2730 return get_result($adjp, $name);
27             }
28 42         148 }
29             }
30              
31             sub initial {
32 20     20 0 55 my $p = shift;
33            
34 20 50 66     69 if(! List::Vectorize::is_array_ref($p) and ! List::Vectorize::is_hash_ref($p)) {
35 0         0 croak "ERROR: P-values should be stored in an array reference or a hash reference.\n";
36             }
37            
38 20         279 my $name = [];
39 20 100       66 if(List::Vectorize::is_hash_ref($p)) {
40 9         64 $name = [ keys %{$p} ];
  9         53  
41 9         23 $p = [ values %{$p} ];
  9         35  
42             }
43            
44 20 50 33     152 if(max($p) > 1 or min($p) < 0) {
45 0         0 croak "ERROR: P-values should between 0 and 1.\n";
46             }
47            
48 20         2779 return ($name, $p);
49             }
50              
51             sub get_result {
52 18     18 0 64 my ($adjp, $name) = @_;
53            
54 18 100       67 if(is_empty($name)) {
55 10         227 return $adjp;
56             }
57             else {
58 8         127 my $result;
59 8         32 for (0..$#$name) {
60 82         234 $result->{$name->[$_]} = $adjp->[$_];
61             }
62 8         97 return $result;
63             }
64             }
65              
66             # R code: pmin(1, n * p)
67             sub _bonferroni {
68 2     2   4 my $p = shift;
69 2         8 my $n = len($p);
70            
71 2         31 my $adjp = multiply($n, $p);
72            
73 2         2935 return pmin(1, $adjp);
74             }
75              
76             # R code: i = 1:n
77             # o = order(p)
78             # ro = order(o)
79             # pmin(1, cummax((n - i + 1) * p[o]))[ro]
80             sub _holm {
81 2     2   6 my $p = shift;
82 2         9 my $n = len($p);
83            
84 2         33 my $i = seq(1, $n);
85 2         188 my $o = order($p);
86 2         281 my $ro = order($o);
87            
88 2         221 my $adjp = multiply(minus($n + 1, $i), subset($p, $o));
89 2         5714 $adjp = cumf($adjp, \&max);
90 2         1021 $adjp = pmin(1, $adjp);
91            
92 2         100 return subset($adjp, $ro);
93             }
94              
95             # R code: i = 1:n
96             # o = order(p)
97             # p = p[o]
98             # ro = order[o]
99             # q = pa = rep.int(min(n * p/i), n)
100             # for (j in (n - 1):2) {
101             # ij = 1:(n - j + 1)
102             # i2 = (n - j + 2):n
103             # q1 <- min(j * p[i2]/(2:j))
104             # q[ij] <- pmin(j * p[ij], q1)
105             # q[i2] <- q[n - j + 1]
106             # pa <- pmax(pa, q)
107             # }
108             # pmax(pa, p)[ro]
109             sub _hommel {
110 2     2   4 my $p = shift;
111 2         7 my $n = len($p);
112            
113 2         26 my $i = seq(1, $n);
114 2         162 my $o = order($p);
115 2         239 $p = subset($p, $o);
116 2         529 my $ro = order($o);
117            
118 2         191 my $pa = rep(min(divide(multiply($n, $p), $i)), $n);
119 2         5622 my $q = copy($pa);
120            
121             # set the first index as 1
122 2         828 unshift(@$p, 0);
123 2         7 unshift(@$q, 0);
124 2         8 unshift(@$pa, 0);
125            
126 2         7 my $ij;
127             my $i2;
128 2         0 my $q1;
129 2         5 for my $j (@{seq($n - 1, 2)}) {
  2         11  
130            
131 16         1051 $ij = seq(1, $n - $j + 1);
132 16         1282 $i2 = seq($n - $j + 2, $n);
133 16         838 $q1 = min(divide(multiply($j, subset($p, $i2)), seq(2, $j)));
134 16         34829 subset_value($q, $ij, pmin(multiply($j, subset($p, $ij)), $q1));
135 16         5181 subset_value($q, $i2, $q->[$n - $j + 1]);
136 16         6061 $pa = pmax($pa, $q);
137             }
138            
139 2         95 shift(@$p);
140 2         6 shift(@$q);
141 2         5 shift(@$pa);
142            
143 2         6 my $adjp = pmax($pa, $p);
144 2         94 return subset($adjp, $ro);
145             }
146              
147             # R code: i = n:1
148             # o <- order(p, decreasing = TRUE)
149             # ro <- order(o)
150             # pmin(1, cummin((n - i + 1) * p[o]))[ro]
151             sub _hochberg {
152            
153 2     2   5 my $p = shift;
154 2         6 my $n = len($p);
155 2         25 my $i = seq($n, 1);
156            
157 2     48   141 my $o = order($p, sub {$_[1] <=> $_[0]});
  48         264  
158 2         44 my $ro = order($o);
159            
160 2         263 my $adjp = multiply(minus($n+1, $i), subset($p, $o));
161 2         5145 $adjp = cumf($adjp, \&min);
162 2         1008 $adjp = pmin(1, $adjp);
163 2         91 return subset($adjp, $ro);
164             }
165              
166             # R code: i <- n:1
167             # o <- order(p, decreasing = TRUE)
168             # ro <- order(o)
169             # pmin(1, cummin(n/i * p[o]))[ro]
170             sub _BH {
171 2     2   8 my $p = shift;
172 2         7 my $n = len($p);
173 2         31 my $i = seq($n, 1);
174            
175 2     48   241 my $o = order($p, sub {$_[1] <=> $_[0]});
  48         286  
176 2         14 my $ro = order($o);
177            
178 2         240 my $adjp = multiply(divide($n, $i), subset($p, $o));
179 2         4858 $adjp = cumf($adjp, \&min);
180 2         1378 $adjp = pmin(1, $adjp);
181 2         105 return subset($adjp, $ro);
182              
183             }
184              
185             # R code: i <- n:1
186             # o <- order(p, decreasing = TRUE)
187             # ro <- order(o)
188             # q <- sum(1/(1L:n))
189             # pmin(1, cummin(q * n/i * p[o]))[ro]
190             sub _BY {
191            
192 2     2   5 my $p = shift;
193 2         6 my $n = len($p);
194 2         24 my $i = seq($n, 1);
195            
196 2     48   138 my $o = order($p, sub {$_[1] <=> $_[0]});
  48         232  
197 2         12 my $ro = order($o);
198            
199 2         187 my $q = sum(divide(1, seq(1, $n)));
200 2         2912 my $adjp = multiply(divide($q*$n, $i), subset($p, $o));
201 2         5520 $adjp = cumf($adjp, \&min);
202 2         1074 $adjp = pmin(1, $adjp);
203 2         135 return subset($adjp, $ro);
204             }
205              
206             sub pmin {
207 26     26 0 21511 my $array1 = shift;
208 26         56 my $array2 = shift;
209            
210 26     188   164 return mapply($array1, $array2, sub {min(\@_)});
  188         38896  
211             }
212              
213             sub pmax {
214 18     18 0 36 my $array1 = shift;
215 18         33 my $array2 = shift;
216            
217 18     196   105 return mapply($array1, $array2, sub {max(\@_)});
  196         29453  
218             }
219              
220              
221             # the default setting can not assure the success of calculation
222             # so this function should run under eval
223             # eval(q(qvalue $p));
224             # if ($@) change_other_settings
225             # @p = qw(0.76 0.30 0.40 0.43 0.27 0.79 0.66 0.36 0.12 0.16 0.52 0.04 0.67 0.44 0.40 0.09 0.51 0.38 0.49 0.68)
226             sub _qvalue {
227 8     8   18 my $p = shift;
228 8         34 my %setup = ('lambda' => multiply(seq(0, 90, 5), 0.01),
229             'robust' => 0,
230             @_ );
231            
232 8         18951 my $lambda = $setup{'lambda'};
233 8         21 my $robust = $setup{'robust'};
234              
235 8 50       30 if(! List::Vectorize::is_array_ref($lambda)) {
236 0         0 croak "ERROR: lambda should be an array reference. If your lambda is a single number, you should set it as an array reference containing one element, such as [0.1].";
237             }
238            
239 8 50 66     78 if(len($lambda) > 1 and len($lambda) < 4) {
240 0         0 croak "ERROR: If length of lambda greater than 1, you need at least 4 values.";
241             }
242 8 50 33     186 if(len($lambda) > 1 and (min($lambda) < 0 or max($lambda) >= 1)) {
      66        
243 0         0 croak "ERROR: Lambda must be within [0, 1).";
244             }
245            
246             # m <- length(p)
247 8         1068 my $m = len($p);
248 8         90 my $pi0;
249            
250 8 100       23 if(len($lambda) == 1) {
251 2         24 $lambda = $lambda->[0];
252 2 50 33     17 if($lambda < 0 or $lambda >= 1) {
253 0         0 croak "ERROR: Lambda must be within [0, 1).";
254             }
255             # pi0 <- mean(p >= lambda)/(1-lambda)
256 2     20   15 $pi0 = mean(test($p, sub {_approximately_compare($_[0], ">=", $lambda)})) / (1 - $lambda);
  20         253  
257             # pi0 <- min(pi0,1)
258 2         315 $pi0 = min(c($pi0, 1));
259             }
260             else {
261             # pi0 <- rep(0,length(lambda))
262 6         72 $pi0 = rep(0, len($lambda));
263             # for(i in 1:length(lambda)) {
264             # pi0[i] <- mean(p >= lambda[i])/(1-lambda[i])
265             # }
266             $pi0 = sapply(seq(0, len($lambda) - 1),
267             sub {
268 114     114   25507 my $current_lambda = $lambda->[$_[0]];
269 114         602 mean(test($p, sub {_approximately_compare($_[0], ">=", $current_lambda)})) / (1 - $current_lambda);
  1197         19747  
270 6         1281 });
271              
272             # minpi0 <- min(pi0)
273 6         1469 my $minpi0 = min($pi0);
274            
275             # mse <- rep(0,length(lambda))
276 6         462 my $mse = rep(0, len($lambda));
277             # pi0.boot <- rep(0,length(lambda))
278 6         1289 my $pi0_boot = rep(0, len($lambda));
279 6         910 for (1..100) {
280             # p.boot <- sample(p,size=m,replace=TRUE)
281 600         999910 my $p_boot = sample($p, $m, 'replace' => 1);
282             # for(i in 1:length(lambda)) {
283             # pi0.boot[i] <- mean(p.boot>lambda[i])/(1-lambda[i])
284             # }
285             $pi0_boot = sapply(seq(0, len($lambda) - 1),
286             sub {
287 11400     11400   2340940 my $current_lambda = $lambda->[$_[0]];
288 11400         49756 mean(test($p_boot, sub {$_[0] > $current_lambda})) / (1 - $current_lambda);
  119700         1816079  
289 600         1432254 });
290             # mse <- mse + (pi0.boot-minpi0)^2
291 600         129036 $mse = plus($mse, power(minus($pi0_boot, $minpi0), 2));
292             }
293            
294             # pi0 <- min(pi0[mse==min(mse)])
295 6         10072 my $min_mse = min($mse);
296            
297 6     114   507 $pi0 = min(subset($pi0, which(test($mse, sub {_approximately_compare($_[0], "==", $min_mse)}))));
  114         1708  
298 6         2331 $pi0 = min(c($pi0, 1));
299             }
300            
301 8 100       1168 if($pi0 <= 0) {
302 2         467 croak "ERROR: The estimated pi0 ($pi0) <= 0. Check that you have valid p-values or use another lambda method.";
303             }
304              
305             # u <- order(p)
306 6         26 my $u = order($p);
307              
308             # this function is not the same to the function in qvalue
309             # but returns the same result
310             # returns number of observations less than or equal to
311             sub qvalue_rank {
312 6     6 0 13 my $x = shift;
313             return sapply($x, sub {
314 62     62   8473 my $current_x = $_[0];
315 62         233 sum(test($x, sub {$_[0] <= $current_x}));
  642         8503  
316 6         29 });
317             }
318            
319             # v <- qvalue.rank(p)
320 6         679 my $v = qvalue_rank($p);
321             # qvalue <- pi0*m*p/v
322 6         989 my $qvalue = divide(multiply($pi0, $m, $p), $v);
323            
324 6 100       19182 if($robust) {
325             # qvalue <- pi0*m*p/(v*(1-(1-p)^m))
326 2         14 $qvalue = divide(multiply($pi0, $m, $p), multiply($v, minus(1, power(minus(1, $p), $m))));
327             }
328            
329             # qvalue[u[m]] <- min(qvalue[u[m]],1)
330 6         6524 $qvalue->[$u->[$#$u]] = min(c($qvalue->[$u->[$#$u]], 1));
331             # for(i in (m-1):1) {
332             # qvalue[u[i]] <- min(qvalue[u[i]],qvalue[u[i+1]],1)
333             # }
334 6         845 for my $i (@{seq($#$u-1, 0)}) {
  6         26  
335 56         2298 $qvalue->[$u->[$i]] = min([$qvalue->[$u->[$i]], $qvalue->[$u->[$i+1]], 1]);
336             }
337            
338 6         235 return $qvalue;
339             }
340              
341             # vectorized power calculation
342             sub power {
343 602     602 0 1370691 my $array = shift;
344 602         1387 my $power = shift;
345            
346 602     11420   3203 my $res = sapply($array, sub {$_[0]**($power)});
  11420         78503  
347 602         6353 return $res;
348             }
349              
350             # approximately
351             sub _approximately_compare {
352 1331     1331   2273 my $left = shift;
353 1331         2226 my $sign = shift;
354 1331         1891 my $right = shift;
355              
356 1331 100 66     4627 if($sign eq ">" or $sign eq ">=") {
    50 33        
    50          
357 1217 100       2918 if($left > $right) {
    100          
358 593         1483 return 1;
359             }
360             elsif(abs($left - $right) < EPS) {
361 38         1420 return 1;
362             }
363             else {
364 586         18941 return 0;
365             }
366             }
367             elsif($sign eq "<" or $sign eq "<=") {
368 0 0       0 if($left < $right) {
    0          
369 0         0 return 1;
370             }
371             elsif(abs($left - $right) < EPS) {
372 0         0 return 1;
373             }
374             else {
375 0         0 return 0;
376             }
377             }
378             elsif($sign eq "==") {
379 114 100       277 if(abs($left - $right) < EPS) {
380 26         1062 return 1;
381             }
382             else {
383 88         2575 return 0;
384             }
385             }
386             }
387              
388             __END__