Skip to content

Commit 6307acb

Browse files
committed
Prim test
1 parent 95a9388 commit 6307acb

File tree

2 files changed

+126
-3
lines changed

2 files changed

+126
-3
lines changed

stan/math/prim/prob/gaussian_copula_cholesky_lpdf.hpp

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -125,15 +125,15 @@ T_return gaussian_copula_cholesky_lpdf(
125125
for (size_t i = 1; i < size_mvt_lcdf_tuple; i++) {
126126
check_size_match(function,
127127
"Size of one of the vectors of "
128-
"the CDF functor tuple",
128+
"the LCDF functor tuple",
129129
math::size(lcdf_tuple_vec[i]),
130130
"Size of the first vector of the "
131131
"location variable",
132132
size_lcdf_tuple);
133133
}
134134

135135
check_size_match(function, "Size of random variable", size_y,
136-
"size of CDF functor tuple", size_lcdf_tuple);
136+
"size of LCDF functor tuple", size_lcdf_tuple);
137137
check_size_match(function, "Size of random variable", size_y,
138138
"rows of covariance parameter", chol.rows());
139139
check_size_match(function, "Size of random variable", size_y,
@@ -149,7 +149,7 @@ T_return gaussian_copula_cholesky_lpdf(
149149
T_return lp(0);
150150
for (size_t i = 0; i < size_vec; i++) {
151151
const auto& u = internal::apply_lcdfs(y_vec[i], lcdf_tuple_vec[i]);
152-
check_bounded(function, "CDF-transformed inputs", u, NEGATIVE_INFTY, 0);
152+
check_bounded(function, "LCDF-transformed inputs", u, NEGATIVE_INFTY, 0);
153153
q[i] = to_vector(std_normal_log_qf(u));
154154
lp -= std_normal_lpdf<propto>(q[i]);
155155
}
Lines changed: 123 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,123 @@
1+
#include <stan/math/prim.hpp>
2+
#include <gtest/gtest.h>
3+
4+
TEST(ProbDistributionsGaussianCopulaCholesky, EqualsMultiNormalCholesky) {
5+
using stan::math::gaussian_copula_cholesky_lpdf;
6+
using stan::math::multi_normal_cholesky_lpdf;
7+
using stan::math::diag_pre_multiply;
8+
9+
Eigen::VectorXd y(2);
10+
y << 1, 3;
11+
12+
Eigen::VectorXd mu(2);
13+
mu << 0.1, 0;
14+
15+
Eigen::VectorXd sd(2);
16+
sd << 2, 1;
17+
18+
Eigen::MatrixXd chol(2, 2);
19+
chol << 1, 0, 0.5, 0.8660254037844386;
20+
21+
auto norm_lcdf = [](auto&& y, auto&& mu, auto&& sigma) {
22+
return stan::math::normal_lcdf(y, mu, sigma);
23+
};
24+
25+
auto std_norm_lcdf = [](auto&& y) {
26+
return stan::math::std_normal_lcdf(y);
27+
};
28+
29+
// y[0] ~ normal(0.1, 2)
30+
// y[1] ~ normal(0, 1)
31+
auto lcdf_functors = std::make_tuple(
32+
std::make_tuple(norm_lcdf, mu[0], sd[0]),
33+
std::make_tuple(std_norm_lcdf)
34+
);
35+
36+
double log_prob =
37+
stan::math::normal_lpdf(y[0], 0.1, 2) +
38+
stan::math::std_normal_lpdf(y[1]) +
39+
gaussian_copula_cholesky_lpdf(y, lcdf_functors, chol);
40+
41+
42+
double expected_log_prob = multi_normal_cholesky_lpdf(
43+
y, mu, diag_pre_multiply(sd, chol));
44+
45+
EXPECT_FLOAT_EQ(log_prob, expected_log_prob);
46+
}
47+
48+
TEST(ProbDistributionsGaussianCopulaCholesky, NonNormalMarginals) {
49+
Eigen::VectorXd y(2);
50+
y << 10, 6;
51+
52+
Eigen::MatrixXd chol(2, 2);
53+
chol << 1, 0, 0.2, 0.9797958971;
54+
55+
auto gamma_lcdf = [](auto&& y, auto&& shape, auto&& scale) {
56+
return stan::math::gamma_lcdf(y, shape, scale);
57+
};
58+
59+
auto exp_lcdf = [](auto&& y, auto&& scale) {
60+
return stan::math::exponential_lcdf(y, scale);
61+
};
62+
63+
// y[0] ~ gamma(2, 1)
64+
// y[1] ~ exponential(2)
65+
auto lcdf_functors = std::make_tuple(
66+
std::make_tuple(gamma_lcdf, 2.0, 1.0),
67+
std::make_tuple(exp_lcdf, 2.0)
68+
);
69+
70+
double log_prob =
71+
stan::math::gamma_lpdf(y[0], 2.0, 1.0) +
72+
stan::math::exponential_lpdf(y[1], 2.0) +
73+
stan::math::gaussian_copula_cholesky_lpdf(y, lcdf_functors, chol);
74+
75+
EXPECT_FLOAT_EQ(log_prob, -16.61005941);
76+
}
77+
78+
TEST(ProbDistributionsGaussianCopulaCholesky, Errors) {
79+
Eigen::VectorXd y(2);
80+
y << 10, 6;
81+
82+
Eigen::VectorXd small_y(1);
83+
small_y << 10;
84+
85+
Eigen::MatrixXd chol(2, 2);
86+
chol << 1, 0, 0.2, 0.9797958971;
87+
88+
auto gamma_lcdf = [](auto&& y, auto&& shape, auto&& scale) {
89+
return stan::math::gamma_lcdf(y, shape, scale);
90+
};
91+
auto exp_lcdf = [](auto&& y, auto&& scale) {
92+
return stan::math::exponential_lcdf(y, scale);
93+
};
94+
95+
auto lcdf_functors = std::make_tuple(
96+
std::make_tuple(gamma_lcdf, 2.0, 1.0),
97+
std::make_tuple(exp_lcdf, 2.0)
98+
);
99+
100+
auto small_lcdf_functors = std::make_tuple(
101+
std::make_tuple(gamma_lcdf, 2.0, 1.0)
102+
);
103+
104+
auto invalid_lcdf_functors = std::make_tuple(
105+
std::make_tuple(gamma_lcdf, 2.0, 1.0),
106+
std::make_tuple([](auto&& y){ return y * 10; })
107+
);
108+
109+
EXPECT_THROW(
110+
stan::math::gaussian_copula_cholesky_lpdf(y, small_lcdf_functors, chol),
111+
std::invalid_argument
112+
);
113+
114+
EXPECT_THROW(
115+
stan::math::gaussian_copula_cholesky_lpdf(small_y, lcdf_functors, chol),
116+
std::invalid_argument
117+
);
118+
119+
EXPECT_THROW(
120+
stan::math::gaussian_copula_cholesky_lpdf(y, invalid_lcdf_functors, chol),
121+
std::domain_error
122+
);
123+
}

0 commit comments

Comments
 (0)