commit be065d5aa2799de2ac4959b868a7fb1b2a654cfd
parent fb2953f658fcf76c9961ed82f1cee58f30689ed6
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Mon, 23 Jul 2018 12:19:52 +0200
Small update of the Microfacet/BSDF/Phase type
Make optional the pdf argument of the sample functions. Make optional
the "type" argument of the BSDF sample function.
Diffstat:
22 files changed, 129 insertions(+), 57 deletions(-)
diff --git a/src/ssf.h b/src/ssf.h
@@ -74,7 +74,7 @@ struct ssf_bsdf_type {
const double wo[3], /* Normalized outgoing direction */
const double N[3], /* Normalized surface normal */
double wi[3], /* Sampled normalized incoming direction */
- int* type, /* Type of the sampled component. Combination of ssf_bxdf_flag */
+ int* type, /* Sampled component. Combination of ssf_bxdf_flag. Can be NULL */
double* pdf); /* PDF to sample wi wrt wo. Can be NULL */
/* Evaluate the BxDF wrt `wo' and `wi' */
@@ -156,7 +156,7 @@ struct ssf_phase_type {
res_T (*init)(struct mem_allocator* allocator, void* bsdf); /*Can be NULL*/
void (*release)(void* bsdf); /* Can be NULL */
- /* Sample a direction `wi' wrt `wo'. TODO comment the PDF */
+ /* Sample a direction `wi' wrt `wo'. */
void
(*sample)
(void* phase,
@@ -172,6 +172,7 @@ struct ssf_phase_type {
const double wo[3], /* Normalized outgoing direction */
const double wi[3]); /* Normalized incoming direction */
+ /* Probability density function */
double
(*pdf)
(void* phase,
@@ -296,14 +297,14 @@ SSF_API const struct ssf_microfacet_distribution_type ssf_blinn_distribution;
SSF_API const struct ssf_microfacet_distribution_type ssf_pillbox_distribution;
/*******************************************************************************
- * Built-in phase functions
+ * Built-in phase functions.
******************************************************************************/
-/* Henyey & Greenstein phase function.
- * TODO add comments */
+/* Henyey & Greenstein phase function normalized to 1; p(wo, wi) == pdf(wo, wi).
+ * p(wo, wi) = 1/(4*PI)* (1-g^2) / (1+g^2-2*g*(-wo.wi))^3/2 */
SSF_API const struct ssf_phase_type ssf_phase_hg;
-/* Rayleigh phase function
- * TODO add comments */
+/* Rayleigh phase function normalized to 1; p(wo, wi) == pdf(wo, wi).
+ * p(wo, wi) = 3/(16*PI) * (1+(-wo.wi)^2) */
SSF_API const struct ssf_phase_type ssf_phase_rayleigh;
/*******************************************************************************
diff --git a/src/ssf_beckmann_distribution.c b/src/ssf_beckmann_distribution.c
@@ -93,7 +93,7 @@ beckmann_distribution_sample
dir[2] = cos_theta;
d33_muld3(wh, d33_basis(basis, N), dir);
ASSERT(d3_dot(wh, N) > 0);
- *pdf = exp(-sin2_theta/(cos2_theta*m2)) / (PI*m2*cos_theta*cos2_theta);
+ if(pdf) *pdf = exp(-sin2_theta/(cos2_theta*m2))/(PI*m2*cos_theta*cos2_theta);
}
static double
diff --git a/src/ssf_blinn_distribution.c b/src/ssf_blinn_distribution.c
@@ -64,7 +64,7 @@ blinn_distribution_sample
double cos_theta, sin_theta;
double phi;
double u, v;
- ASSERT(distrib && rng && N && wh && pdf);
+ ASSERT(distrib && rng && N && wh);
ASSERT(d3_is_normalized(N));
u = ssp_rng_canonical(rng);
@@ -78,8 +78,9 @@ blinn_distribution_sample
dir[2] = cos_theta;
d33_muld3(dir, d33_basis(basis, N), dir);
- *pdf = (blinn->exponent+2)/(2*PI)*pow(cos_theta, blinn->exponent+1);
d3_set(wh, dir);
+
+ if(pdf) *pdf = (blinn->exponent+2)/(2*PI)*pow(cos_theta, blinn->exponent+1);
}
static double
diff --git a/src/ssf_bsdf.c b/src/ssf_bsdf.c
@@ -125,11 +125,10 @@ ssf_bsdf_sample
int* type,
double* out_pdf)
{
- double R, pdf;
- ASSERT(bsdf && rng && wo && N && wi && type);
+ double R;
+ ASSERT(bsdf && rng && wo && N && wi);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N));
- R = bsdf->type.sample(bsdf->data, rng, wo, N, wi, type, &pdf);
- if(out_pdf) *out_pdf = pdf;
+ R = bsdf->type.sample(bsdf->data, rng, wo, N, wi, type, out_pdf);
return R;
}
diff --git a/src/ssf_lambertian_reflection.c b/src/ssf_lambertian_reflection.c
@@ -51,13 +51,13 @@ lambertian_reflection_sample
double* pdf)
{
double sample[3];
- ASSERT(data && rng && wo && N && wi && type && pdf);
+ ASSERT(data && rng && wo && N && wi);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N) && d3_dot(wo, N) > 0);
(void)wo;
ssp_ran_hemisphere_cos(rng, N, sample, pdf);
d3_set(wi, sample);
- *type = SSF_REFLECTION | SSF_DIFFUSE;
+ if(type) *type = SSF_REFLECTION | SSF_DIFFUSE;
return ((struct lambertian_reflection*)data)->reflectivity;
}
diff --git a/src/ssf_microfacet_distribution.c b/src/ssf_microfacet_distribution.c
@@ -128,11 +128,9 @@ ssf_microfacet_distribution_sample
double wh[3],
double* out_pdf)
{
- double pdf;
ASSERT(distrib && rng && N && wh);
ASSERT(d3_is_normalized(N));
- distrib->type.sample(distrib->data, rng, N, wh, &pdf);
- if(out_pdf) *out_pdf = pdf;
+ distrib->type.sample(distrib->data, rng, N, wh, out_pdf);
}
double
diff --git a/src/ssf_microfacet_reflection.c b/src/ssf_microfacet_reflection.c
@@ -83,14 +83,14 @@ microfacet_reflection_sample
double wh[3];
double pdf_wh;
double R;
- ASSERT(data && wo && N && wi && pdf && type);
+ ASSERT(data && wo && N && wi);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N) && d3_dot(wo, N) > 0);
ASSERT(bsdf->distrib && bsdf->fresnel);
ssf_microfacet_distribution_sample(bsdf->distrib, rng, N, wh, &pdf_wh);
reflect(dir, wo, wh);
- *pdf = pdf_wh / (4.0*fabs(d3_dot(wo, N)));
- *type = SSF_REFLECTION | SSF_GLOSSY;
+ if(pdf) *pdf = pdf_wh / (4.0*fabs(d3_dot(wo, N)));
+ if(type) *type = SSF_REFLECTION | SSF_GLOSSY;
R = d3_dot(dir, N) > 0 ? ssf_fresnel_eval(bsdf->fresnel, d3_dot(dir, wh)) : 0;
d3_set(wi, dir);
return R;
@@ -105,7 +105,7 @@ microfacet_reflection_pdf
double pdf_wh;
ASSERT(data && wo && N && wi);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N) && d3_is_normalized(wi));
- ASSERT(d3_dot(wo, N) > 0 && d3_dot(wi, N) > 0);
+ ASSERT(d3_dot(wo, N) > 0);
ASSERT(bsdf->distrib);
d3_normalize(wh, d3_add(wh, wi, wo));
if(d3_dot(wh, N) < 0) d3_minus(wh, wh);
@@ -131,7 +131,7 @@ microfacet2_reflection_sample
double wh[3];
double p;
double troughput = 1;
- ASSERT(data && wo && N && wi && pdf && type);
+ ASSERT(data && wo && N && wi);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N) && d3_dot(N, wo) > 0.0);
ASSERT(bsdf->distrib);
@@ -155,8 +155,8 @@ microfacet2_reflection_sample
reflect(dir, dir, wh);
}
- *pdf = NaN;
- *type = SSF_REFLECTION | SSF_GLOSSY;
+ if(pdf) *pdf = NaN;
+ if(type) *type = SSF_REFLECTION | SSF_GLOSSY;
d3_set(wi, dir);
return troughput;
}
diff --git a/src/ssf_phase.c b/src/ssf_phase.c
@@ -123,11 +123,9 @@ ssf_phase_sample
double wi[3],
double* out_pdf)
{
- double pdf;
ASSERT(phase && rng && wo && wi);
ASSERT(d3_is_normalized(wo));
- phase->type.sample(phase->data, rng, wo, wi, &pdf);
- if(out_pdf) *out_pdf = pdf;
+ phase->type.sample(phase->data, rng, wo, wi, out_pdf);
}
double
diff --git a/src/ssf_pillbox_distribution.c b/src/ssf_pillbox_distribution.c
@@ -82,7 +82,8 @@ pillbox_distribution_sample
dir[1] = sin(phi) * sin_theta;
dir[2] = cos_theta;
d33_muld3(wh, d33_basis(basis, N), dir);
- *pdf = cos_theta / (PI-PI*cos2_theta_max);
+
+ if(pdf) *pdf = cos_theta / (PI-PI*cos2_theta_max);
}
static double
diff --git a/src/ssf_specular_dielectric_dielectric_interface.c b/src/ssf_specular_dielectric_dielectric_interface.c
@@ -63,7 +63,7 @@ ssf_specular_dielectric_dielectric_interface_sample
double cos_wo_N;
double eta; /* Ratio of eta_i / eta_t */
double R;
- ASSERT(bsdf && rng && wi && N && wo && type && pdf);
+ ASSERT(bsdf && rng && wi && N && wo);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N));
ASSERT(d3_dot(wo, N) > -1.e-6);
(void)rng;
@@ -82,15 +82,15 @@ ssf_specular_dielectric_dielectric_interface_sample
SSF(fresnel_dielectric_dielectric_setup(fresnel, bsdf->eta_i, bsdf->eta_t));
R = ssf_fresnel_eval(fresnel, cos_wo_N);
- *pdf = INF;
+ if(pdf) *pdf = INF;
/* Sample the output direction wrt R */
if(ssp_rng_canonical(rng) < R) {
reflect(wi, wo, N);
- *type = SSF_SPECULAR | SSF_REFLECTION;
+ if(type) *type = SSF_SPECULAR | SSF_REFLECTION;
} else {
d3_set(wi, refracted);
- *type = SSF_SPECULAR | SSF_TRANSMISSION;
+ if(type) *type = SSF_SPECULAR | SSF_TRANSMISSION;
}
return 1;
}
diff --git a/src/ssf_specular_reflection.c b/src/ssf_specular_reflection.c
@@ -47,14 +47,14 @@ specular_reflection_sample
struct specular_reflection* brdf = data;
double cos_wo_N;
double cos_wi_N;
- ASSERT(rng && wi && N && wo && type && pdf);
+ ASSERT(rng && wi && N && wo);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N) && d3_dot(wo, N) > 0);
(void)rng;
/* Reflect the outgoing direction wo with respect to the surface normal N */
reflect(wi, wo, N);
- *pdf = INF;
- *type = SSF_REFLECTION | SSF_SPECULAR;
+ if(pdf) *pdf = INF;
+ if(type) *type = SSF_REFLECTION | SSF_SPECULAR;
cos_wi_N = cos_wo_N = d3_dot(wo, N);
return ssf_fresnel_eval(brdf->fresnel, cos_wi_N);
diff --git a/src/ssf_thin_specular_dielectric.c b/src/ssf_thin_specular_dielectric.c
@@ -66,7 +66,7 @@ thin_specular_dielectric_sample
double R, T;
double rho1, rho2, rho2_sqr;
double tau, tau_sqr;
- ASSERT(bsdf && rng && wi && N && wo && type && pdf);
+ ASSERT(bsdf && rng && wi && N && wo);
ASSERT(d3_is_normalized(wo) && d3_is_normalized(N));
ASSERT(d3_dot(wo, N) > -1.e-6);
(void)rng;
@@ -74,8 +74,8 @@ thin_specular_dielectric_sample
eta = bsdf->eta_i / bsdf->eta_t;
if(!refract(wt, wo, N, eta)) { /* Total reflection */
reflect(wi, wo, N);
- *pdf = INF;
- *type = SSF_SPECULAR | SSF_REFLECTION;
+ if(pdf) *pdf = INF;
+ if(type) *type = SSF_SPECULAR | SSF_REFLECTION;
return 1;
}
fresnel = bsdf->fresnel;
@@ -103,15 +103,15 @@ thin_specular_dielectric_sample
}
#endif
- *pdf = INF;
+ if(pdf) *pdf = INF;
/* Importance sample the BTDF wrt R */
if(ssp_rng_uniform_double(rng, 0, R + T) < R) {
reflect(wi, wo, N);
- *type = SSF_SPECULAR | SSF_REFLECTION;
+ if(type) *type = SSF_SPECULAR | SSF_REFLECTION;
} else { /* Sample the transmissive part */
d3_minus(wi, wo);
- *type = SSF_SPECULAR | SSF_TRANSMISSION;
+ if(type) *type = SSF_SPECULAR | SSF_TRANSMISSION;
}
if(bsdf->absorption == 0) {
/* avoid numerical loss of energy if no absorption */
diff --git a/src/test_ssf_bsdf.c b/src/test_ssf_bsdf.c
@@ -67,8 +67,8 @@ bsdf_sample
CHK(d3_eq(BxDF->wo, wo) == 1);
CHK(d3_eq(BxDF->N, N) == 1);
d3(wi, 1.0, 2.0, 3.0);
- *type = 314;
- *pdf = 4;
+ if(type) *type = 314;
+ if(pdf) *pdf = 4;
return BxDF->reflectivity;
}
@@ -190,10 +190,12 @@ main(int argc, char** argv)
data->rng = rng;
data->reflectivity = 0.1234;
+ i = 0, pdf = 0;
CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, &i, NULL) == 0.1234);
- CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, &i, &pdf) == 0.1234);
CHK(i == 314);
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, &pdf) == 0.1234);
CHK(pdf == 4);
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, NULL) == 0.1234);
data->reflectivity = 0.314;
CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, &i, &pdf) == 0.314);
diff --git a/src/test_ssf_lambertian_reflection.c b/src/test_ssf_lambertian_reflection.c
@@ -66,6 +66,13 @@ main(int argc, char** argv)
CHK(eq_eps(pdf, d3_dot(wi, N)/PI, 1.e-6) == 1);
CHK(type == (SSF_DIFFUSE|SSF_REFLECTION));
+ pdf = 0, type = 0;
+ CHK(ssf_bsdf_sample(brdf, rng, wo, N, wi, NULL, &pdf) == R);
+ CHK(eq_eps(pdf, d3_dot(wi, N)/PI, 1.e-6) == 1);
+ CHK(ssf_bsdf_sample(brdf, rng, wo, N, wi, &type, NULL) == R);
+ CHK(type == (SSF_DIFFUSE|SSF_REFLECTION));
+ CHK(ssf_bsdf_sample(brdf, rng, wo, N, wi, NULL, NULL) == R);
+
d3_normalize(wo, d3(wo, 1.0, 0.0, 1.0));
N[0] = ssp_rng_uniform_double(rng, -1, 1);
N[1] = ssp_rng_uniform_double(rng, -1, 1);
diff --git a/src/test_ssf_microfacet_distribution.c b/src/test_ssf_microfacet_distribution.c
@@ -60,11 +60,10 @@ ufacet_sample
CHK(ufacet != NULL);
CHK(N != NULL);
CHK(wh != NULL);
- CHK(pdf != NULL);
CHK(ufacet->rng == rng);
CHK(d3_eq(ufacet->N, N) == 1);
d3_normalize(wh, d3(wh, 1.0, 2.0, 3.0));
- *pdf = ufacet->pdf;
+ if(pdf) *pdf = ufacet->pdf;
}
static double
diff --git a/src/test_ssf_microfacet_reflection.c b/src/test_ssf_microfacet_reflection.c
@@ -30,7 +30,11 @@ main(int argc, char** argv)
struct ssp_rng* rng;
double N[3];
double wo[3];
- size_t i, NSTEPS=100000;
+ double wi[3];
+ double pdf;
+ double R;
+ size_t i, NSTEPS=10000;
+ int type;
(void)argc, (void)argv;
mem_init_proxy_allocator(&allocator, &mem_default_allocator);
@@ -56,11 +60,46 @@ main(int argc, char** argv)
ssp_ran_sphere_uniform(rng, N, NULL);
ssp_ran_hemisphere_cos(rng, N, wo, NULL);
FOR_EACH(i, 0, NSTEPS) {
- double wi[3], pdf;
- int type;
CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, &type, &pdf) == 1);
CHK(type == (SSF_GLOSSY | SSF_REFLECTION));
+ CHK(IS_NaN(pdf));
CHK(d3_dot(wi, N) > 0);
+
+ pdf = 0, type = 0;
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, &pdf) == 1);
+ CHK(IS_NaN(pdf));
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, &type, NULL) == 1);
+ CHK(type == (SSF_GLOSSY | SSF_REFLECTION));
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, NULL) == 1);
+ }
+
+ CHK(ssf_bsdf_ref_put(bsdf) == RES_OK);
+ CHK(ssf_bsdf_create(&allocator, &ssf_microfacet_reflection, &bsdf) == RES_OK);
+
+ CHK(ssf_microfacet_reflection_setup(NULL, NULL, NULL) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(bsdf, NULL, NULL) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(NULL, F, NULL) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(bsdf, F, NULL) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(NULL, NULL, D) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(bsdf, NULL, D) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(NULL, F, D) == RES_BAD_ARG);
+ CHK(ssf_microfacet_reflection_setup(bsdf, F, D) == RES_OK);
+
+ FOR_EACH(i, 0, NSTEPS) {
+ R = ssf_bsdf_sample(bsdf, rng, wo, N, wi, &type, &pdf);
+ CHK(R == (d3_dot(wi, N) > 0 ? 1 : 0));
+ CHK(type == (SSF_GLOSSY | SSF_REFLECTION));
+ CHK(eq_eps(pdf, ssf_bsdf_pdf(bsdf, wo, N, wi), 1.e-6));
+
+ pdf = 0, type = 0;
+ R = ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, &pdf);
+ CHK(R == (d3_dot(wi, N) > 0 ? 1 : 0));
+ CHK(eq_eps(pdf, ssf_bsdf_pdf(bsdf, wo, N, wi), 1.e-6));
+ R = ssf_bsdf_sample(bsdf, rng, wo, N, wi, &type, NULL);
+ CHK(R == (d3_dot(wi, N) > 0 ? 1 : 0));
+ CHK(type == (SSF_GLOSSY | SSF_REFLECTION));
+ R = ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, NULL);
+ CHK(R == (d3_dot(wi, N) > 0 ? 1 : 0));
}
CHK(ssf_bsdf_ref_put(bsdf) == RES_OK);
diff --git a/src/test_ssf_phase.c b/src/test_ssf_phase.c
@@ -59,7 +59,7 @@ phase_sample
CHK(phase->rng == rng);
CHK(d3_eq(phase->wo, wo) == 1);
d3(wi, 1.0, 2.0, 3.0);
- *pdf = phase->pdf;
+ if(pdf) *pdf = phase->pdf;
}
static double
@@ -166,6 +166,7 @@ main(int argc, char** argv)
data->value = 0.1234;
data->pdf = 4;
+ pdf = 0;
ssf_phase_sample(phase, rng, data->wo, wi, NULL);
ssf_phase_sample(phase, rng, data->wo, wi, &pdf);
CHK(pdf == 4);
diff --git a/src/test_ssf_phase_hg.c b/src/test_ssf_phase_hg.c
@@ -26,6 +26,8 @@ main(int argc, char** argv)
struct ssf_phase* phase;
struct ssf_phase* dummy;
struct ssp_rng* rng;
+ double wo[3];
+ double wi[3];
size_t i;
(void)argc, (void)argv;
@@ -45,7 +47,6 @@ main(int argc, char** argv)
FOR_EACH(i, 0, NG) {
const double g = ssp_rng_uniform_double(rng, -1, +1);
- double wo[3];
double sum_cos = 0;
double sum_cos_sqr = 0;
double E, V, SE;
@@ -55,7 +56,6 @@ main(int argc, char** argv)
ssp_ran_sphere_uniform(rng, wo, NULL);
FOR_EACH(isamp, 0, NSAMPS) {
- double wi[3];
double w[3];
double weight;
double pdf;
@@ -82,6 +82,10 @@ main(int argc, char** argv)
CHK(eq_eps(E, g, 3*SE));
}
+ ssp_ran_sphere_uniform(rng, wo, NULL);
+ ssf_phase_sample(phase, rng, wo, wi, NULL);
+ CHK(d3_is_normalized(wi));
+
CHK(ssf_phase_ref_put(phase) == RES_OK);
CHK(ssf_phase_ref_put(dummy) == RES_OK);
CHK(ssp_rng_ref_put(rng) == RES_OK);
diff --git a/src/test_ssf_phase_rayleigh.c b/src/test_ssf_phase_rayleigh.c
@@ -25,6 +25,7 @@ main(int argc, char** argv)
struct ssf_phase* phase;
struct ssp_rng* rng;
double wo[3];
+ double wi[3];
double sum_cos = 0;
double sum_cos_sqr = 0;
double sum_cos_sqr2 = 0;
@@ -39,7 +40,6 @@ main(int argc, char** argv)
ssp_ran_sphere_uniform(rng, wo, NULL);
FOR_EACH(i, 0, NSAMPS) {
- double wi[3];
double w[3];
double weight;
double pdf;
@@ -64,13 +64,16 @@ main(int argc, char** argv)
E = sum_cos / NSAMPS;
V = sum_cos_sqr / NSAMPS - E*E;
SE = sqrt(V/NSAMPS);
- CHK(eq_eps(E, 0, SE));
+ CHK(eq_eps(E, 0, SE*3));
/* On average cos(-wo,wi)^2 should be 2/5 */
E = sum_cos_sqr / NSAMPS;
V = sum_cos_sqr2 / NSAMPS - E*E;
SE = sqrt(V/NSAMPS);
- CHK(eq_eps(E, 2.0/5.0, SE));
+ CHK(eq_eps(E, 2.0/5.0, SE*3));
+
+ ssf_phase_sample(phase, rng, wo, wi, NULL);
+ CHK(d3_is_normalized(wi));
CHK(ssf_phase_ref_put(phase) == RES_OK);
CHK(ssp_rng_ref_put(rng) == RES_OK);
diff --git a/src/test_ssf_specular_reflection.c b/src/test_ssf_specular_reflection.c
@@ -56,6 +56,13 @@ main(int argc, char** argv)
CHK(IS_INF(pdf) == 1);
CHK(type == (SSF_SPECULAR|SSF_REFLECTION));
+ pdf = 0, type = 0;
+ CHK(ssf_bsdf_sample(brdf, rng, wo, N, wi, NULL, &pdf) == R);
+ CHK(IS_INF(pdf));
+ CHK(ssf_bsdf_sample(brdf, rng, wo, N, wi, &type, NULL) == R);
+ CHK(type == (SSF_SPECULAR|SSF_REFLECTION));
+ CHK(ssf_bsdf_sample(brdf, rng, wo, N, wi, NULL, NULL) == R);
+
d3_normalize(wo, d3(wo, 1.0, 1.0, 0.0));
CHK(ssf_specular_reflection_setup(brdf, fresnel) == RES_OK);
R = ssf_bsdf_sample(brdf, rng, wo, N, wi, &type, &pdf);
diff --git a/src/test_ssf_thin_specular_dielectric.c b/src/test_ssf_thin_specular_dielectric.c
@@ -115,6 +115,15 @@ main(int argc, char** argv)
CHK(IS_INF(pdf) == 1);
CHK(d3_eq_eps(wi, d3(tmp, -wo[0], wo[1], 0), 1.e-6) == 1);
CHK(type == (SSF_SPECULAR | SSF_REFLECTION));
+
+ /* Check optional arguments */
+ pdf = 0, type = 0;
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, &pdf) == R);
+ CHK(IS_INF(pdf));
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, &type, NULL) == R);
+ CHK(type == (SSF_SPECULAR | SSF_REFLECTION));
+ CHK(ssf_bsdf_sample(bsdf, rng, wo, N, wi, NULL, NULL) == R);
+
FOR_EACH(i, 0, NSTEPS) {
R = ssf_bsdf_sample(bsdf, rng, wo, N, wi, &type, &pdf);
CHK(type & SSF_REFLECTION);
diff --git a/src/test_ssf_utils.h b/src/test_ssf_utils.h
@@ -282,6 +282,9 @@ check_microfacet_distribution
D = ssf_microfacet_distribution_eval(distrib, N, wh);
weight = D * d3_dot(wh, N) / pdf;
CHK(eq_eps(weight, 1.0, 1.e-6) == 1);
+ CHK(eq_eps(pdf, ssf_microfacet_distribution_pdf(distrib, N, wh), 1.e-6f));
+
+ ssf_microfacet_distribution_sample(distrib, rng, N, wh, NULL);
}
}