9.8 Measured BSDFs

The reflection models introduced up to this point represent index of refraction changes at smooth and rough boundaries, which constitute the basic building blocks of surface appearance. More complex materials (e.g., paint on primer, metal under a layer of enamel) can sometimes be approximated using multiple interfaces with participating media between them; the layered material model presented in Section 14.3 is based on that approach.

However, many real-world materials are beyond the scope of even such layered models. Examples include:

  • Materials characterized by wave-optical phenomena that produce striking directionally varying coloration. Examples include iridescent paints, insect wings, and holographic paper.
  • Materials with rough interfaces. In pbrt, we have chosen to model such surfaces using microfacet theory and the Trowbridge–Reitz distribution. However, it is important to remember that both of these are models that generally do not match real-world behavior perfectly.
  • Surfaces with non-standard microstructure. For example, a woven fabric composed of two different yarns looks like a surface from a distance, but its directional intensity and color variation are not well-described by any standard BRDF model due to the distinct reflectance properties of fiber-based microgeometry.

Instead of developing numerous additional specialized BxDFs, we will now pursue another way of reproducing such challenging materials in a renderer: by interpolating measurements of real-world material samples to create a data-driven reflectance model. The resulting MeasuredBxDF only models surface reflection, though the approach can in principle generalize to transmission as well.

<<MeasuredBxDF Definition>>= 
class MeasuredBxDF { public: <<MeasuredBxDF Public Methods>> 
PBRT_CPU_GPU MeasuredBxDF(const MeasuredBxDFData *brdf, const SampledWavelengths &lambda) : brdf(brdf), lambda(lambda) {} static MeasuredBxDFData *BRDFDataFromFile(const std::string &filename, Allocator alloc); PBRT_CPU_GPU SampledSpectrum f(Vector3f wo, Vector3f wi, TransportMode mode) const; PBRT_CPU_GPU pstd::optional<BSDFSample> Sample_f(Vector3f wo, Float uc, Point2f u, TransportMode mode, BxDFReflTransFlags sampleFlags) const; PBRT_CPU_GPU Float PDF(Vector3f wo, Vector3f wi, TransportMode mode, BxDFReflTransFlags sampleFlags) const; PBRT_CPU_GPU void Regularize() {} PBRT_CPU_GPU static constexpr const char *Name() { return "MeasuredBxDF"; } std::string ToString() const; PBRT_CPU_GPU BxDFFlags Flags() const { return (BxDFFlags::Reflection | BxDFFlags::Glossy); }
private: <<MeasuredBxDF Private Methods>> 
static Float theta2u(Float theta) { return std::sqrt(theta * (2 / Pi)); } static Float phi2u(Float phi) { return phi * (1 / (2 * Pi)) + .5f; } static Float u2theta(Float u) { return Sqr(u) * (Pi / 2.f); } static Float u2phi(Float u) { return (2.f * u - 1.f) * Pi; }
<<MeasuredBxDF Private Members>> 
const MeasuredBxDFData *brdf; SampledWavelengths lambda;
};

Measuring reflection in a way that is practical while producing information in a form that is convenient for rendering is a challenging problem. We begin by explaining these challenges for motivation.

Consider the task of measuring the BRDF of a sheet of brushed aluminum: we could illuminate a sample of the material from a set of n incident directions left-parenthesis theta Subscript normal i Superscript left-parenthesis k right-parenthesis Baseline comma phi Subscript normal i Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis with k equals 1 comma ellipsis comma n and use some kind of sensor (e.g., a photodiode) to record the reflected light scattered along a set of m outgoing directions left-parenthesis theta Subscript normal o Superscript left-parenthesis l right-parenthesis Baseline comma phi Subscript normal o Superscript left-parenthesis l right-parenthesis Baseline right-parenthesis with l equals 1 comma ellipsis comma m . These n times m measurements could be stored on disk and interpolated at runtime to approximate the BRDF at intermediate configurations left-parenthesis theta Subscript normal i Baseline comma phi Subscript normal i Baseline comma theta Subscript normal o Baseline comma phi Subscript normal o Baseline right-parenthesis . However, closer scrutiny of such an approach reveals several problems:

Figure 9.33: Specialized Hardware for BSDF Acquisition. The term goniophotometer (or gonioreflectometer) refers to a typically motorized platform that can simultaneously illuminate and observe a material sample from arbitrary pairs of directions. The device on the left (at Cornell University, built by Cyberware Inc., image courtesy of Steve Marschner) rotates camera (2 degrees of freedom) and light arms (1 degree of freedom) around a centered sample pedestal that can also rotate about its vertical axis. The device on the right (at EPFL, built by pab advanced technologies Ltd) instead uses a static light source and a rotating sensor arm (2 degrees of freedom). The vertical material sample holder then provides 2 rotational degrees of freedom to cover the full 4D domain of the BSDF.

  • BSDFs of polished materials are highly directionally peaked. Perturbing the incident or outgoing direction by as little as 1 degree can change the measured reflectance by orders of magnitude. This implies that the set of incident and outgoing directions must be sampled fairly densely.
  • Accurate positioning in spherical coordinates is difficult to perform by hand and generally requires mechanical aids. For this reason, such measurements are normally performed using a motorized gantry known as a goniophotometer or gonioreflectometer. Figure 9.33 shows two examples of such machines. Light stages consisting of a rigid assembly of hundreds of LEDs around a sample are sometimes used to accelerate measurement, though at the cost of reduced directional resolution.
  • Sampling each direction using a 1 degree spacing in spherical coordinates requires roughly one billion sample points. Storing gigabytes of measurement data is possible but undesirable, yet the time that would be spent for a full measurement is even more problematic: assuming that the goniophotometer can reach a configuration left-parenthesis theta Subscript normal i Baseline comma phi Subscript normal i Baseline comma theta Subscript normal o Baseline comma phi Subscript normal o Baseline right-parenthesis within 1 second (a reasonable estimate for the devices shown in Figure 9.33), over 34 years of sustained operation would be needed to measure a single material.

In sum, the combination of high-frequency signals, the 4D domain of the BRDF, and the curse of dimensionality conspire to make standard measurement approaches all but impractical.

Figure 9.34: Adaptive BRDF Sample Placement. (a) Regular sampling of the incident and outgoing directions is a prohibitively expensive way of measuring and storing BRDFs due to the curse of dimensionality (here, only 2 of the 4 dimensions are shown). (b) A smaller number of samples can yield a more accurate interpolant if their placement is informed by the material’s reflectance behavior.

While there is no general antidote against the curse of dimensionality, the previous example involving a dense discretization of the 4D domain is clearly excessive. For example, peaked BSDFs that concentrate most of their energy into a small set of angles tend to be relatively smooth away from the peak. Figure 9.34 shows how a more specialized sample placement that is informed by the principles of specular reflection can drastically reduce the number of sample points that are needed to obtain a desired level of accuracy. Figure 9.35 shows how the roughness of the surface affects the desired distribution of samples—for example, smooth surfaces allow sparse sampling outside of the specular lobe.

Figure 9.35: The Effect of Surface Roughness on Adaptive BRDF Sample Placement. The two plots visualize BRDF values of two materials with different roughnesses for varying directions omega Subscript normal i and fixed omega Subscript normal o . Circles indicate adaptively chosen measurement locations, which are used to create the interpolant implemented in the MeasuredBxDF class. (a) The measurement locations broadly cover the hemisphere given a relatively rough material. (b) For a more specular material the samples are concentrated in the region around the specular peak. Changing the outgoing direction moves the specular peak; hence the sample locations must depend on omega Subscript normal o .

The MeasuredBxDF therefore builds on microfacet theory and the distribution of visible normals to create a more efficient physically informed sampling pattern. The rationale underlying this choice is that while microfacet theory may not perfectly predict the reflectance of a specific material, it can at least approximately represent how energy is (re-)distributed throughout the 4D domain. Applying it enables the use of a relatively coarse set of measurement locations that suffice to capture the function’s behavior.

Concretely, the method works by transforming regular grid points using visible normal sampling (Section 9.6.4) and performing a measurement at each sampled position. If the microfacet sampling routine is given by a function upper R colon script upper S squared times left-bracket 0 comma 1 right-bracket squared right-arrow script upper S squared and u Superscript left-parenthesis k right-parenthesis with k equals 1 comma ellipsis comma n denotes input samples arranged on a grid covering the 2D unit square, then we have a sequence of measurements script upper M Superscript left-parenthesis k right-parenthesis :

script upper M Superscript left-parenthesis k right-parenthesis Baseline equals f Subscript normal r Baseline left-parenthesis omega Subscript normal o Baseline comma upper R left-parenthesis omega Subscript normal o Baseline comma u Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis right-parenthesis comma

where f Subscript normal r Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline right-parenthesis refers to the real-world BRDF of a material sample, as measured by a goniophotometer (or similar device) in directions omega Subscript normal o and omega Subscript normal i Baseline equals upper R left-parenthesis omega Subscript normal o Baseline comma u Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis . This process must be repeated for different values of omega Subscript normal o to also capture variation in the other direction argument. Evaluating the BRDF requires the inverse upper R Superscript negative 1 of the transformation, which yields a position on left-bracket 0 comma 1 right-bracket squared that can be used to interpolate the measurements script upper M Superscript left-parenthesis k right-parenthesis . Figure 9.36 illustrates both directions of this mapping.

Figure 9.36: Visible Normal Sampling as a Parameterization. The MeasuredBxDF leverages visible normal sampling as a parameterization upper R of the unit sphere. Here, it is used in a deterministic fashion to transform a set of grid points u Superscript left-parenthesis k right-parenthesis with k equals 1 comma ellipsis comma n into spherical directions to be measured. Evaluating the resulting BRDF representation requires the inverse upper R Superscript negative 1 followed by linear interpolation within the regular grid of measurements.

This procedure raises several questions: first, the non-random use of a method designed for Monte Carlo sampling may be unexpected. To see why this works, remember that the inversion method (Section 2.3) evaluates the inverse of a distribution’s cumulative distribution function (CDF). Besides being convenient for sampling, this inverse CDF can also be interpreted as a parameterization of the target domain from the unit square. This parameterization smoothly warps the domain so that regions with a high contribution occupy a correspondingly large amount of the unit square. The MeasuredBxDF then simply measures and stores BRDF values in these “improved” coordinates. Note that the material does not have to agree with microfacet theory for this warping to be valid, though the sampling pattern is much less efficient and requires a denser discretization when the material’s behavior deviates significantly.

Another challenge is that parameterization guiding the measurement requires a microfacet approximation of the material, but such an approximation would normally be derived from an existing measurement. We will shortly show how to resolve this chicken-and-egg problem and assume for now that a suitable model is available.

Measurement Through a Parameterization

A flaw of the reparameterized measurement sequence in Equation (9.41) is that the values script upper M Superscript left-parenthesis k right-parenthesis may differ by many orders of magnitude, which means that simple linear interpolation between neighboring data points is unlikely to give satisfactory results. We instead use the following representation that transforms measurements prior to storage in an effort to reduce this variation:

script upper M Superscript left-parenthesis k right-parenthesis Baseline equals StartFraction f Subscript normal r Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis cosine theta Subscript normal i Superscript left-parenthesis k right-parenthesis Baseline Over p left-parenthesis omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis EndFraction comma

where omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis Baseline equals upper R left-parenthesis omega Subscript normal o Baseline comma u Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis , and p left-parenthesis omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis Baseline right-parenthesis denotes the density of direction omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis according to visible normal sampling.

If f Subscript normal r was an analytic BRDF (e.g., a microfacet model) and u Superscript left-parenthesis k right-parenthesis a 2D uniform variate, then Equation (9.42) would simply be the weight of a Monte Carlo importance sampling strategy, typically with a carefully designed mapping upper R and density p that make this weight near-constant to reduce variance.

In the present context, f Subscript normal r represents real-world data, while p and upper R encapsulate a microfacet approximation of the material under consideration. We therefore expect script upper M Superscript left-parenthesis k right-parenthesis to take on near-constant values when the material is well-described by a microfacet model, and more marked deviations otherwise. This can roughly be interpreted as measuring the difference (in a multiplicative sense) between the real world and the microfacet simplification. Figure 9.37 visualizes the effect of the transformation in Equation (9.42).

Figure 9.37: Reparameterized BRDF Visualization. This figure illustrates the representation of two material samples: a metallic sample swatch from the L3-37 robot in the film Solo: A Star Wars Story (Walt Disney Studios Motion Pictures) and a pearlescent vehicle vinyl wrap (TeckWrap International Inc.). Each column represents a measurement of a separate outgoing direction omega Subscript normal o . For both materials, the first row visualizes the measured directions omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis . The subsequent row plots the “raw” reparameterized BRDF of Equation (9.41), where each pixel represents one of the grid points u Superscript left-parenthesis k right-parenthesis Baseline element-of left-bracket 0 comma 1 right-bracket squared identified with omega Subscript normal i Baseline Superscript left-parenthesis k right-parenthesis . The final row shows transformed measurements corresponding to Equation (9.42) that are more uniform and easier to interpolate. Note that these samples are both isotropic, which is why a few measurements for different elevation angles suffice. In the anisotropic case, the left-parenthesis theta Subscript normal o Baseline comma phi Subscript normal o Baseline right-parenthesis domain must be covered more densely.

BRDF Evaluation

Evaluating the data-driven BRDF requires the inverse of these steps. Suppose that script upper M left-parenthesis dot right-parenthesis implements an interpolation based on the grid of measurement points script upper M Superscript left-parenthesis k right-parenthesis . Furthermore, suppose that we have access to the inverse upper R Superscript negative 1 Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline right-parenthesis that returns the “random numbers” u that would cause visible normal sampling to generate a particular incident direction (i.e., upper R left-parenthesis omega Subscript normal o Baseline comma u right-parenthesis equals omega Subscript normal i ). Accessing script upper M left-parenthesis dot right-parenthesis through upper R Superscript negative 1 then provides a spherical interpolation of the measurement data.

We must additionally multiply by the density p left-parenthesis omega Subscript normal i Baseline right-parenthesis , and divide by the cosine factor to undo corresponding transformations introduced in Equation (9.42), which yields the final form of the data-driven BRDF:

f Subscript normal r Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline right-parenthesis equals StartFraction script upper M left-parenthesis upper R Superscript negative 1 Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline right-parenthesis right-parenthesis p left-parenthesis omega Subscript normal i Baseline right-parenthesis Over cosine theta Subscript normal i Baseline EndFraction period

Generalized Microfacet Model

A major difference between the microfacet model underlying the ConductorBxDF and the approximation used here is that we replace the Trowbridge–Reitz model with an arbitrary data-driven microfacet distribution. This improves the model’s ability to approximate the material being measured. At the same time, it implies that previously used simplifications and analytic solutions are no longer available and must be revisited.

We begin with the Torrance–Sparrow sampling density from Equation (9.28),

p left-parenthesis omega Subscript normal i Baseline right-parenthesis equals StartFraction upper D Subscript omega Sub Subscript normal o Subscript Baseline left-parenthesis omega Subscript normal m Baseline right-parenthesis Over 4 left-parenthesis omega Subscript normal o Baseline dot omega Subscript normal m Baseline right-parenthesis EndFraction comma

which references the visible normal sampling density upper D Subscript omega Sub Subscript Baseline left-parenthesis omega Subscript normal m Baseline right-parenthesis from Equation (9.23). Substituting the definition of the masking function from Equation (9.18) into upper D Subscript omega Sub Subscript normal o Baseline left-parenthesis omega Subscript normal m Baseline right-parenthesis and rearranging terms yields

upper D Subscript omega Sub Subscript normal o Subscript Baseline left-parenthesis omega Subscript normal m Baseline right-parenthesis equals StartFraction 1 Over sigma left-parenthesis omega Subscript normal o Baseline right-parenthesis EndFraction upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis max left-parenthesis 0 comma omega Subscript normal o Baseline dot omega Subscript normal m Baseline right-parenthesis comma

where

sigma left-parenthesis omega Subscript normal o Baseline right-parenthesis equals integral Underscript script upper H squared left-parenthesis bold n Subscript Baseline right-parenthesis Endscripts upper D left-parenthesis omega Subscript Baseline right-parenthesis max left-parenthesis 0 comma omega Subscript normal o Baseline dot omega Subscript Baseline right-parenthesis normal d omega Subscript Baseline

provides a direction-dependent normalization of the visible normal distribution. For valid reflection directions ( omega Subscript normal o Baseline dot omega Subscript normal m Baseline greater-than 0 ), the PDF of generated samples then simplifies to

p left-parenthesis omega Subscript normal i Baseline right-parenthesis equals StartFraction upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis Over 4 sigma left-parenthesis omega Subscript normal o Baseline right-parenthesis EndFraction period

Substituting this density into the BRDF from Equation (9.43) produces

f Subscript normal r Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline right-parenthesis equals StartFraction script upper M left-parenthesis upper R Superscript negative 1 Baseline left-parenthesis omega Subscript normal o Baseline comma omega Subscript normal i Baseline right-parenthesis right-parenthesis upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis Over 4 sigma left-parenthesis omega Subscript normal o Baseline right-parenthesis cosine theta Subscript normal i Baseline EndFraction period

The MeasuredBxDF implements this expression using data-driven representations of upper D left-parenthesis dot right-parenthesis and sigma left-parenthesis dot right-parenthesis .

Finding the Initial Microfacet Model

We finally revisit the chicken-and-egg problem from before: practical measurement using the presented approach requires a suitable microfacet model—specifically, a microfacet distribution upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis . Yet it remains unclear how this distribution could be obtained without access to an already existing BRDF measurement.

The key idea to resolve this conundrum is that the microfacet distribution upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis is a 2D quantity, which means that it remains mostly unaffected by the curse of dimensionality. Acquiring this function is therefore substantially cheaper than a measurement of the full 4D BRDF.

Suppose that the material being measured perfectly matches microfacet theory in the sense that it is described by the Torrance–Sparrow BRDF from Equation (9.33). Then we can measure the material’s retroreflection (i.e., omega Subscript normal i Baseline equals omega Subscript normal o Baseline equals omega Subscript ), which is given by

f Subscript normal r Baseline left-parenthesis normal p Subscript Baseline comma omega Subscript Baseline comma omega Subscript Baseline right-parenthesis equals StartFraction upper D left-parenthesis omega Subscript Baseline right-parenthesis upper F left-parenthesis omega Subscript Baseline dot omega Subscript Baseline right-parenthesis upper G left-parenthesis omega Subscript Baseline comma omega Subscript Baseline right-parenthesis Over 4 cosine squared theta EndFraction proportional-to StartFraction upper D left-parenthesis omega Subscript Baseline right-parenthesis upper G 1 left-parenthesis omega Subscript Baseline right-parenthesis Over cosine squared theta EndFraction period

The last step of the above equation removes constant terms including the Fresnel reflectance and introduces the reasonable assumption that shadowing/masking is perfectly correlated given omega Subscript normal i Baseline equals omega Subscript normal o and thus occurs only once. Substituting the definition of upper G 1 from Equation (9.18) and rearranging yields the following relationship of proportionality:

upper D left-parenthesis omega Subscript Baseline right-parenthesis proportional-to f Subscript normal r Baseline left-parenthesis normal p Subscript Baseline comma omega Subscript Baseline comma omega Subscript Baseline right-parenthesis cosine theta integral Underscript script upper H squared left-parenthesis bold n Subscript Baseline right-parenthesis Endscripts upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis max left-parenthesis 0 comma omega Subscript Baseline dot omega Subscript normal m Baseline right-parenthesis normal d omega Subscript normal m Baseline period

This integral equation can be solved by measuring f Subscript normal r Baseline left-parenthesis normal p Subscript Baseline comma omega Subscript Baseline Subscript j Baseline comma omega Subscript Baseline Subscript j Baseline right-parenthesis for n directions omega Subscript Baseline Subscript j and using those measurements for initial guesses of upper D left-parenthesis omega Subscript Baseline Subscript j Baseline right-parenthesis . A more accurate estimate of upper D can then be found using an iterative solution procedure where the estimated values of upper D are used to estimate the integrals on the right hand side of Equation (9.46) for all of the omega Subscript Baseline Subscript j s. This process quickly converges within a few iterations.

9.8.1 Basic Data Structures

MeasuredBxDFData holds data pertaining to reflectance measurements and the underlying parameterization. Because the data for an isotropic BRDF is typically a few megabytes and the data for an anisotropic BRDF may be over 100, each measured BRDF that is used in the scene is stored in memory only once. As instances of MeasuredBxDF are created at surface intersections during rendering, they can then store just a pointer to the appropriate MeasuredBxDFData. Code not included here adds the ability to initialize instances of this type from binary .bsdf files containing existing measurements.

<<MeasuredBxDFData Definition>>= 
struct MeasuredBxDFData { <<MeasuredBxDFData Public Members>> 
pstd::vector<float> wavelengths; PiecewiseLinear2D<3> spectra; PiecewiseLinear2D<0> ndf; PiecewiseLinear2D<2> vndf; PiecewiseLinear2D<0> sigma; bool isotropic; PiecewiseLinear2D<2> luminance; MeasuredBxDFData(Allocator alloc) : ndf(alloc), sigma(alloc), vndf(alloc), luminance(alloc), spectra(alloc), wavelengths(alloc) {} static MeasuredBxDFData *Create(const std::string &filename, Allocator alloc); std::string ToString() const { return StringPrintf("[ MeasuredBxDFData filename: %s ]", filename); } std::string filename;
};

Measured BRDFs are represented by spectral measurements at a set of discrete wavelengths that are stored in wavelengths. The actual measurements are stored in spectra.

<<MeasuredBxDFData Public Members>>= 
pstd::vector<float> wavelengths; PiecewiseLinear2D<3> spectra;

The template class PiecewiseLinear2D represents a piecewise-linear interpolant on the 2D unit square with samples arranged on a regular grid. The details of its implementation are relatively technical and reminiscent of other interpolants in this book; hence we only provide an overview of its capabilities and do not include its full implementation here.

The class is parameterized by a Dimension template parameter that extends the 2D interpolant to higher dimensions—for example, PiecewiseLinear2D<1> stores a 3D grid of values, and PiecewiseLinear2D<3> used above for spectra is a 5D quantity. The class provides three key methods:

template <size_t Dimension> class PiecewiseLinear2D { public: Float Evaluate(Point2f pos, Float... params); PLSample Sample(Point2f u, Float... params); PLSample Invert(Point2f p, Float... params); };

where PLSample is defined as

struct PLSample { Point2f p; Float pdf; };

Evaluate() takes a position monospace p monospace o monospace s element-of left-bracket 0 comma 1 right-bracket squared and then additional Float parameters to perform a lookup using multidimensional linear interpolation according to the value of Dimension.

Sample() warps u element-of left-bracket 0 comma 1 right-bracket squared via inverse transform sampling (i.e., proportional to the stored linear interpolant), returning both the sampled position on left-bracket 0 comma 1 right-bracket squared and the associated density as a PLSample. The additional parameters passed via params are used as conditional variables that restrict sampling to a 2D slice of a higher-dimensional function. For example, invoking the method PiecewiseLinear2D<3>::Sample() with a uniform 2D variate and parameters 0.1 , 0.2 , and 0.3 would importance sample the 2D slice upper I left-parenthesis 0.1 comma 0.2 comma 0.3 comma dot comma dot right-parenthesis of a pentalinear interpolant upper I .

Finally, Invert() implements the exact inverse of Sample(). Invoking it with the position computed by Sample() will recover the input u value up to rounding error.

Additional PiecewiseLinear1D instances are used to (redundantly) store the normal distribution upper D left-parenthesis omega Subscript normal m Baseline right-parenthesis in ndf, the visible normal distribution upper D Subscript omega Sub Subscript normal o Baseline left-parenthesis omega Subscript normal m Baseline right-parenthesis parameterized by omega Subscript normal o Baseline equals left-parenthesis theta Subscript normal o Baseline comma phi Subscript normal o Baseline right-parenthesis in vndf, and the normalization constant sigma left-parenthesis omega Subscript normal o Baseline right-parenthesis in sigma. The data structure also records whether the material is isotropic, in which case the dimensionality of some of the piecewise-linear interpolants can be reduced.

<<MeasuredBxDFData Public Members>>+=  
PiecewiseLinear2D<0> ndf; PiecewiseLinear2D<2> vndf; PiecewiseLinear2D<0> sigma; bool isotropic;

9.8.2 Evaluation

Figure 9.38: Three Measured BRDFs. (a) TeckWrap Amber Citrine vinyl wrapping film. (b) Purple acrylic felt. (c) Silk from an Indian sari with two colors of yarn.

Following these preliminaries, we can now turn to evaluating the measured BRDF for a pair of directions. See Figure 9.38 for examples of the variety of types of reflection that the measured representation can reproduce.

The only information that must be stored as MeasuredBxDF member variables in order to implement the BxDF interface methods is a pointer to the BRDF measurement data and the set of wavelengths at which the BRDF is to be evaluated.

<<MeasuredBxDF Private Members>>= 
const MeasuredBxDFData *brdf; SampledWavelengths lambda;

BRDF evaluation then follows the approach described in Equation (9.45).

<<MeasuredBxDF Method Definitions>>= 
SampledSpectrum MeasuredBxDF::f(Vector3f wo, Vector3f wi, TransportMode mode) const { <<Check for valid reflection configurations>> 
if (!SameHemisphere(wo, wi)) return SampledSpectrum(0); if (wo.z < 0) { wo = -wo; wi = -wi; }
<<Determine half-direction vector omega Subscript normal m >> 
Vector3f wm = wi + wo; if (LengthSquared(wm) == 0) return SampledSpectrum(0); wm = Normalize(wm);
<<Map omega Subscript normal o and omega Subscript normal m to the unit square left-bracket 0 comma 1 right-bracket squared >> 
Float theta_o = SphericalTheta(wo), phi_o = std::atan2(wo.y, wo.x); Float theta_m = SphericalTheta(wm), phi_m = std::atan2(wm.y, wm.x); Point2f u_wo(theta2u(theta_o), phi2u(phi_o)); Point2f u_wm(theta2u(theta_m), phi2u(brdf->isotropic ? (phi_m - phi_o) : phi_m)); u_wm[1] = u_wm[1] - pstd::floor(u_wm[1]);
<<Evaluate inverse parameterization upper R Superscript negative 1 >> 
PLSample ui = brdf->vndf.Invert(u_wm, phi_o, theta_o);
<<Evaluate spectral 5D interpolant>> 
SampledSpectrum fr; for (int i = 0; i < NSpectrumSamples; ++i) fr[i] = std::max<Float>(0, brdf->spectra.Evaluate(ui.p, phi_o, theta_o, lambda[i]));
<<Return measured BRDF value>> 
return fr * brdf->ndf.Evaluate(u_wm) / (4 * brdf->sigma.Evaluate(u_wo) * CosTheta(wi));
}

Zero reflection is returned if the specified directions correspond to transmission through the surface. Otherwise, the directions omega Subscript normal i and omega Subscript normal o are mirrored onto the positive hemisphere if necessary.

<<Check for valid reflection configurations>>= 
if (!SameHemisphere(wo, wi)) return SampledSpectrum(0); if (wo.z < 0) { wo = -wo; wi = -wi; }

The next code fragment determines the associated microfacet normal and handles an edge case that occurs in near-grazing configurations.

<<Determine half-direction vector omega Subscript normal m >>= 
Vector3f wm = wi + wo; if (LengthSquared(wm) == 0) return SampledSpectrum(0); wm = Normalize(wm);

A later step requires that omega Subscript normal o and omega Subscript normal m are mapped onto the unit square left-bracket 0 comma 1 right-bracket squared , which we do in two steps: first, by converting the directions to spherical coordinates, which are then further transformed by helper methods theta2u() and phi2u().

In the isotropic case, the mapping used for omega Subscript normal m subtracts phi Subscript normal o from phi Subscript normal m , which allows the stored tables to be invariant to rotation about the surface normal. This may cause the second dimension of u_wm to fall out of the left-bracket 0 comma 1 right-bracket interval; a subsequent correction fixes this using the periodicity of the azimuth parameter.

<<Map omega Subscript normal o and omega Subscript normal m to the unit square left-bracket 0 comma 1 right-bracket squared >>= 
Float theta_o = SphericalTheta(wo), phi_o = std::atan2(wo.y, wo.x); Float theta_m = SphericalTheta(wm), phi_m = std::atan2(wm.y, wm.x); Point2f u_wo(theta2u(theta_o), phi2u(phi_o)); Point2f u_wm(theta2u(theta_m), phi2u(brdf->isotropic ? (phi_m - phi_o) : phi_m)); u_wm[1] = u_wm[1] - pstd::floor(u_wm[1]);

The two helper functions encapsulate an implementation detail of the storage representation. The function phi2u() uniformly maps left-bracket negative pi comma pi right-bracket onto left-bracket 0 comma 1 right-bracket , while theta2u() uses a nonlinear transformation that places more resolution near theta almost-equals 0 to facilitate storing the microfacet distribution of highly specular materials.

<<MeasuredBxDF Private Methods>>= 
static Float theta2u(Float theta) { return std::sqrt(theta * (2 / Pi)); } static Float phi2u(Float phi) { return phi * (1 / (2 * Pi)) + .5f; }

With this information at hand, we can now evaluate the inverse parameterization to determine the sample values ui.p that would cause visible normal sampling to generate the current incident direction (i.e., upper R left-parenthesis omega Subscript normal o Baseline comma u right-parenthesis equals omega Subscript normal i ).

<<Evaluate inverse parameterization upper R Superscript negative 1 >>= 
PLSample ui = brdf->vndf.Invert(u_wm, phi_o, theta_o);

This position is then used to evaluate a 5D linear interpolant parameterized by the fractional 2D position monospace u monospace i monospace period monospace p element-of left-bracket 0 comma 1 right-bracket squared on the reparameterized incident hemisphere, phi Subscript normal o , theta Subscript normal o , and the wavelength in nanometers. The interpolant must be evaluated once per sample of SampledSpectrum.

<<Evaluate spectral 5D interpolant>>= 
SampledSpectrum fr; for (int i = 0; i < NSpectrumSamples; ++i) fr[i] = std::max<Float>(0, brdf->spectra.Evaluate(ui.p, phi_o, theta_o, lambda[i]));

Finally, fr must be scaled to undo the transformations that were applied to the data to improve the quality of the interpolation and to reduce the required measurement density, giving the computation that corresponds to Equation (9.45).

<<Return measured BRDF value>>= 
return fr * brdf->ndf.Evaluate(u_wm) / (4 * brdf->sigma.Evaluate(u_wo) * CosTheta(wi));

In principle, implementing the Sample_f() and PDF() methods required by the BxDF interface is straightforward: the Sample_f() method could evaluate the forward mapping upper R to perform visible normal sampling based on the measured microfacet distribution using PiecewiseLinear2D::Sample(), and PDF() could evaluate the associated sampling density from Equation (9.44). However, a flaw of such a basic sampling scheme is that the transformed BRDF measurements from Equation (9.42) are generally nonuniform on left-bracket 0 comma 1 right-bracket squared , which can inject unwanted variance into the rendering process. The implementation therefore uses yet another reparameterization based on a luminance tensor that stores the product integral of the spectral dimension of MeasuredBxDFData::spectra and the CIE Y color matching curve.

<<MeasuredBxDFData Public Members>>+= 
PiecewiseLinear2D<2> luminance;

The actual BRDF sampling routine then consists of two steps. First it converts a uniformly distributed sample on left-bracket 0 comma 1 right-bracket squared into another sample u element-of left-bracket 0 comma 1 right-bracket squared that is distributed according to the luminance of the reparameterized BRDF. Following this, visible normal sampling transforms u into a sampled direction omega Subscript normal i and a sampling weight that will have near-constant luminance. Apart from this step, the implementations of Sample_f() and PDF() are similar to the evaluation method and therefore not included here.