added MicrofacetDistribution::computeTransmissionProbability

metadata
Wenzel Jakob 2011-07-11 16:38:46 +02:00
parent 1debcf3c0b
commit 0803cba093
2 changed files with 78 additions and 14 deletions

View File

@ -421,19 +421,21 @@ public:
}
/**
* \brief Precompute the aggregate Fresnel transmittance through
* a rough interface
* \brief Compute a spline representation for the overall Fresnel
* transmittance through a rough interface
*
* This function efficiently computes the integral of
* \int_{S^2} D(m) * (1 - Fr(m<->wi)) dm
* This function essentially computes the integral of
* 1 - \int_{S^2} f(w_i, w_o) * dw_o
* for incident directions 'wi' with a range of different inclinations
* (where Fr denotes the fresnel reflectance). It returns a cubic spline
* interpolation parameterized by the cosine of the angle between 'wi'
* and the (macro-) surface normal.
* (where f denotes a Cook-Torrance style reflectance model). It returns
* a cubic spline interpolation parameterized by the cosine of the angle
* between 'wi' and the (macro-) surface normal.
*
* \remark This only works for isotropic microfacet distributions
*/
CubicSpline *getRoughTransmittance(Float extIOR, Float intIOR, Float alpha, size_t resolution) const {
CubicSpline *computeRoughTransmittance(Float extIOR, Float intIOR, Float alpha, size_t resolution) const {
if (isAnisotropic())
SLog(EError, "MicrofacetDistribution::getRoughTransmission(): only "
SLog(EError, "MicrofacetDistribution::computeRoughTransmission(): only "
"supports isotropic distributions!");
NDIntegrator integrator(1, 2, 5000, 0, 1e-5f);
@ -449,7 +451,7 @@ public:
integral = 0, error = 0;
integrator.integrateVectorized(
boost::bind(&MicrofacetDistribution::integrand, this,
boost::bind(&MicrofacetDistribution::integrand1, this,
wi, extIOR, intIOR, alpha, _1, _2, _3),
min, max, &integral, &error, &nEvals
);
@ -465,6 +467,52 @@ public:
return spline;
}
/**
* \brief Compute a spline representation that gives the probability
* of sampling a transmission event e.g. in the plugin 'roughdielectric'.
*
* Like \ref computeRoughTransmittance, the spline is parameterized by the
* cosine of the angle between the indident direction and the (macro-)
* surface normal.
*
* \remark This function only works for isotropic microfacet distributions
*/
CubicSpline *computeTransmissionProbability(Float extIOR, Float intIOR, Float alpha, size_t resolution) const {
if (isAnisotropic())
SLog(EError, "MicrofacetDistribution::computeTransmissionProbability(): only "
"supports isotropic distributions!");
NDIntegrator integrator(1, 2, 5000, 0, 1e-5f);
CubicSpline *spline = new CubicSpline(resolution);
size_t nEvals, nEvalsTotal = 0;
ref<Timer> timer = new Timer();
Float stepSize = (1.0f-2*Epsilon)/(resolution-1);
for (size_t i=0; i<resolution; ++i) {
Float z = stepSize * i + Epsilon;
Vector wi(std::sqrt(std::max((Float) 0, 1-z*z)), 0, z);
Float min[2] = {0, 0}, max[2] = {1, 1},
integral = 0, error = 0;
integrator.integrateVectorized(
boost::bind(&MicrofacetDistribution::integrand2, this,
wi, extIOR, intIOR, alpha, _1, _2, _3),
min, max, &integral, &error, &nEvals
);
spline->append(z, integral);
nEvalsTotal += nEvals;
}
SLog(EInfo, "Created a " SIZE_T_FMT "-node cubic spline approximation to the "
"transmission probability (integration took %i ms and " SIZE_T_FMT
" function evaluations)", resolution, timer->getMilliseconds(),
nEvalsTotal);
spline->build();
return spline;
}
std::string toString() const {
switch (m_type) {
case EBeckmann: return "beckmann"; break;
@ -477,13 +525,13 @@ public:
}
}
protected:
/// Integrand helper function called by \ref getRoughTransmission
void integrand(const Vector &wi, Float extIOR, Float intIOR, Float alpha,
/// Integrand helper function called by \ref computeRoughTransmission
void integrand1(const Vector &wi, Float extIOR, Float intIOR, Float alpha,
size_t nPts, const Float *in, Float *out) const {
for (int i=0; i<(int) nPts; ++i) {
Normal m = sample(Point2(in[2*i], in[2*i+1]), alpha);
Vector wo = 2 * dot(wi, m) * Vector(m) - wi;
if (Frame::cosTheta(wi) <= 0 || Frame::cosTheta(wo) <= 0) {
if (Frame::cosTheta(wo) <= 0) {
out[i] = 0;
continue;
}
@ -495,6 +543,22 @@ protected:
}
}
/// Integrand helper function called by \ref computeTransmissionProbability
void integrand2(const Vector &wi, Float extIOR, Float intIOR, Float alpha,
size_t nPts, const Float *in, Float *out) const {
for (int i=0; i<(int) nPts; ++i) {
Normal m = sample(Point2(in[2*i], in[2*i+1]), alpha);
Vector wo = 2 * dot(wi, m) * Vector(m) - wi;
if (Frame::cosTheta(wo) <= 0) {
out[i] = 0;
continue;
}
/* Calculate the specular reflection component */
out[i] = 1 - fresnel(dot(wi, m), extIOR, intIOR);
}
}
protected:
EType m_type;
};

View File

@ -130,7 +130,7 @@ public:
if (m_roughTransmittance == NULL) {
Float alpha = m_distribution.transformRoughness(m_alpha->getValue(Intersection()).average());
m_roughTransmittance = m_distribution.getRoughTransmittance(m_extIOR, m_intIOR, alpha, 200);
m_roughTransmittance = m_distribution.computeRoughTransmittance(m_extIOR, m_intIOR, alpha, 200);
}
BSDF::configure();