1 // NOTE: Throughout, we use the symbol ⊙ for convolution.
2 // Since all of our signals are symmetrical, discrete correlation and convolution
3 // is the same operation, and so we won't make a difference in notation.
9 #include <Eigen/Cholesky>
11 #include "deconvolution_sharpen_effect.h"
15 using namespace Eigen;
17 DeconvolutionSharpenEffect::DeconvolutionSharpenEffect()
20 gaussian_radius(0.0f),
24 last_circle_radius(-1.0f),
25 last_gaussian_radius(-1.0f),
26 last_correlation(-1.0f),
29 register_int("matrix_size", &R);
30 register_float("circle_radius", &circle_radius);
31 register_float("gaussian_radius", &gaussian_radius);
32 register_float("correlation", &correlation);
33 register_float("noise", &noise);
36 std::string DeconvolutionSharpenEffect::output_fragment_shader()
39 sprintf(buf, "#define R %u\n", R);
42 assert(R <= 25); // Same limit as Refocus.
45 return buf + read_file("deconvolution_sharpen_effect.frag");
50 // Integral of sqrt(r² - x²) dx over x=0..a.
51 float circle_integral(float a, float r)
58 return 0.25f * M_PI * r * r;
60 return 0.5f * (a * sqrt(r*r - a*a) + r*r * asin(a / r));
63 // Yields the impulse response of a circular blur with radius r.
64 // We basically look at each element as a square centered around (x,y),
65 // and figure out how much of its area is covered by the circle.
66 float circle_impulse_response(int x, int y, float r)
69 // Degenerate case: radius = 0 yields the impulse response.
70 return (x == 0 && y == 0) ? 1.0f : 0.0f;
73 // Find the extents of this cell. Due to symmetry, we can cheat a bit
74 // and pretend we're always in the upper-right quadrant, except when
75 // we're right at an axis crossing (x = 0 or y = 0), in which case we
76 // simply use the evenness of the function; shrink the cell, make
77 // the calculation, and down below we'll normalize by the cell's area.
78 float min_x, max_x, min_y, max_y;
83 min_x = abs(x) - 0.5f;
84 max_x = abs(x) + 0.5f;
90 min_y = abs(y) - 0.5f;
91 max_y = abs(y) + 0.5f;
93 assert(min_x >= 0.0f && max_x >= 0.0f);
94 assert(min_y >= 0.0f && max_y >= 0.0f);
96 float cell_height = max_y - min_y;
97 float cell_width = max_x - min_x;
99 if (min_x * min_x + min_y * min_y > r * r) {
100 // Lower-left corner is outside the circle, so the entire cell is.
103 if (max_x * max_x + max_y * max_y < r * r) {
104 // Upper-right corner is inside the circle, so the entire cell is.
108 // OK, so now we know the cell is partially covered by the circle:
118 // The edge of the circle is defined by x² + y² = r²,
119 // or x = sqrt(r² - y²) (since x is nonnegative).
120 // Find out where the curve crosses our given y values.
121 float mid_x1 = (max_y >= r) ? min_x : sqrt(r * r - max_y * max_y);
122 float mid_x2 = sqrt(r * r - min_y * min_y);
123 if (mid_x1 < min_x) {
126 if (mid_x2 > max_x) {
129 assert(mid_x1 >= min_x);
130 assert(mid_x2 >= mid_x1);
131 assert(max_x >= mid_x2);
133 // The area marked A in the figure above.
134 float covered_area = cell_height * (mid_x1 - min_x);
136 // The area marked B in the figure above. Note that the integral gives the entire
137 // shaded space down to zero, so we need to subtract the rectangle that does not
138 // belong to our cell.
139 covered_area += circle_integral(mid_x2, r) - circle_integral(mid_x1, r);
140 covered_area -= min_y * (mid_x2 - mid_x1);
142 assert(covered_area <= cell_width * cell_height);
143 return covered_area / (cell_width * cell_height);
146 // Compute a ⊙ b. Note that we compute the “full” convolution,
147 // ie., our matrix will be big enough to hold every nonzero element of the result.
148 MatrixXf convolve(const MatrixXf &a, const MatrixXf &b)
150 MatrixXf result(a.rows() + b.rows() - 1, a.cols() + b.cols() - 1);
151 for (int yr = 0; yr < result.rows(); ++yr) {
152 for (int xr = 0; xr < result.cols(); ++xr) {
155 // Given that x_b = x_r - x_a, find the values of x_a where
156 // x_a is in [0, a_cols> and x_b is in [0, b_cols>. (y is similar.)
158 // The second demand gives:
160 // 0 <= x_r - x_a < b_cols
161 // 0 >= x_a - x_r > -b_cols
162 // x_r >= x_a > x_r - b_cols
163 int ya_min = yr - b.rows() + 1;
165 int xa_min = xr - b.rows() + 1;
168 // Now fit to the first demand.
169 ya_min = std::max<int>(ya_min, 0);
170 ya_max = std::min<int>(ya_max, a.rows() - 1);
171 xa_min = std::max<int>(xa_min, 0);
172 xa_max = std::min<int>(xa_max, a.cols() - 1);
174 assert(ya_max >= ya_min);
175 assert(xa_max >= xa_min);
177 for (int ya = ya_min; ya <= ya_max; ++ya) {
178 for (int xa = xa_min; xa <= xa_max; ++xa) {
179 sum += a(ya, xa) * b(yr - ya, xr - xa);
183 result(yr, xr) = sum;
189 // Similar to convolve(), but instead of assuming every element outside
190 // of b is zero, we make no such assumption and instead return only the
191 // elements where we know the right answer. (This is the only difference
193 // This is the same as conv2(a, b, 'valid') in Octave.
195 // a must be the larger matrix of the two.
196 MatrixXf central_convolve(const MatrixXf &a, const MatrixXf &b)
198 assert(a.rows() >= b.rows());
199 assert(a.cols() >= b.cols());
200 MatrixXf result(a.rows() - b.rows() + 1, a.cols() - b.cols() + 1);
201 for (int yr = b.rows() - 1; yr < result.rows() + b.rows() - 1; ++yr) {
202 for (int xr = b.cols() - 1; xr < result.cols() + b.cols() - 1; ++xr) {
205 // Given that x_b = x_r - x_a, find the values of x_a where
206 // x_a is in [0, a_cols> and x_b is in [0, b_cols>. (y is similar.)
208 // The second demand gives:
210 // 0 <= x_r - x_a < b_cols
211 // 0 >= x_a - x_r > -b_cols
212 // x_r >= x_a > x_r - b_cols
213 int ya_min = yr - b.rows() + 1;
215 int xa_min = xr - b.rows() + 1;
218 // Now fit to the first demand.
219 ya_min = std::max<int>(ya_min, 0);
220 ya_max = std::min<int>(ya_max, a.rows() - 1);
221 xa_min = std::max<int>(xa_min, 0);
222 xa_max = std::min<int>(xa_max, a.cols() - 1);
224 assert(ya_max >= ya_min);
225 assert(xa_max >= xa_min);
227 for (int ya = ya_min; ya <= ya_max; ++ya) {
228 for (int xa = xa_min; xa <= xa_max; ++xa) {
229 sum += a(ya, xa) * b(yr - ya, xr - xa);
233 result(yr - b.rows() + 1, xr - b.cols() + 1) = sum;
239 void print_matrix(const MatrixXf &m)
241 for (int y = 0; y < m.rows(); ++y) {
242 for (int x = 0; x < m.cols(); ++x) {
243 printf("%7.4f ", m(x, y));
251 void DeconvolutionSharpenEffect::update_deconvolution_kernel()
253 // Figure out the impulse response for the circular part of the blur.
254 MatrixXf circ_h(2 * R + 1, 2 * R + 1);
255 for (int y = -R; y <= R; ++y) {
256 for (int x = -R; x <= R; ++x) {
257 circ_h(y + R, x + R) = circle_impulse_response(x, y, circle_radius);
261 // Same, for the Gaussian part of the blur. We make this a lot larger
262 // since we're going to convolve with it soon, and it has infinite support
263 // (see comments for central_convolve()).
264 MatrixXf gaussian_h(4 * R + 1, 4 * R + 1);
265 for (int y = -2 * R; y <= 2 * R; ++y) {
266 for (int x = -2 * R; x <= 2 * R; ++x) {
268 if (gaussian_radius < 1e-3) {
269 val = (x == 0 && y == 0) ? 1.0f : 0.0f;
271 val = exp(-(x*x + y*y) / (2.0 * gaussian_radius * gaussian_radius));
273 gaussian_h(y + 2 * R, x + 2 * R) = val;
277 // h, the (assumed) impulse response that we're trying to invert.
278 MatrixXf h = central_convolve(gaussian_h, circ_h);
279 assert(h.rows() == 2 * R + 1);
280 assert(h.cols() == 2 * R + 1);
282 // Normalize the impulse response.
284 for (int y = 0; y < 2 * R + 1; ++y) {
285 for (int x = 0; x < 2 * R + 1; ++x) {
289 for (int y = 0; y < 2 * R + 1; ++y) {
290 for (int x = 0; x < 2 * R + 1; ++x) {
295 // r_uu, the (estimated/assumed) autocorrelation of the input signal (u).
296 // The signal is modelled a standard autoregressive process with the
297 // given correlation coefficient.
299 // We have to take a bit of care with the size of this matrix.
300 // The pow() function naturally has an infinite support (except for the
301 // degenerate case of correlation=0), but we have to chop it off
302 // somewhere. Since we convolve it with a 4*R+1 large matrix below,
303 // we need to make it twice as big as that, so that we have enough
304 // data to make r_vv valid. (central_convolve() effectively enforces
305 // that we get at least the right size.)
306 MatrixXf r_uu(8 * R + 1, 8 * R + 1);
307 for (int y = -4 * R; y <= 4 * R; ++y) {
308 for (int x = -4 * R; x <= 4 * R; ++x) {
309 r_uu(x + 4 * R, y + 4 * R) = pow(correlation, hypot(x, y));
313 // Estimate r_vv, the autocorrelation of the output signal v.
314 // Since we know that v = h ⊙ u and both are symmetrical,
315 // convolution and correlation are the same, and
316 // r_vv = v ⊙ v = (h ⊙ u) ⊙ (h ⊙ u) = (h ⊙ h) ⊙ r_uu.
317 MatrixXf r_vv = central_convolve(r_uu, convolve(h, h));
318 assert(r_vv.rows() == 4 * R + 1);
319 assert(r_vv.cols() == 4 * R + 1);
321 // Similarly, r_uv = u ⊙ v = u ⊙ (h ⊙ u) = h ⊙ r_uu.
322 MatrixXf r_uu_center = r_uu.block(2 * R, 2 * R, 4 * R + 1, 4 * R + 1);
323 MatrixXf r_uv = central_convolve(r_uu_center, h);
324 assert(r_uv.rows() == 2 * R + 1);
325 assert(r_uv.cols() == 2 * R + 1);
327 // Add the noise term (we assume the noise is uncorrelated,
328 // so it only affects the central element).
329 r_vv(2 * R, 2 * R) += noise;
331 // Now solve the Wiener-Hopf equations to find the deconvolution kernel g.
332 // Most texts show this only for the simpler 1D case:
334 // [ r_vv(0) r_vv(1) r_vv(2) ... ] [ g(0) ] [ r_uv(0) ]
335 // [ r_vv(-1) r_vv(0) ... ] [ g(1) ] = [ r_uv(1) ]
336 // [ r_vv(-2) ... ] [ g(2) ] [ r_uv(2) ]
337 // [ ... ] [ g(3) ] [ r_uv(3) ]
339 // (Since r_vv is symmetrical, we can drop the minus signs.)
341 // Generally, row i of the matrix contains (dropping _vv for brevity):
343 // [ r(0-i) r(1-i) r(2-i) ... ]
345 // However, we have the 2D case. We flatten the vectors out to
346 // 1D quantities; this means we must think of the row number
347 // as a pair instead of as a scalar. Row (i,j) then contains:
349 // [ r(0-i,0-j) r(1-i,0-j) r(2-i,0-j) ... r(0-i,1-j) r_(1-i,1-j) r(2-i,1-j) ... ]
351 // g and r_uv are flattened in the same fashion.
353 // Note that even though this matrix is block Toeplitz, it is _not_ Toeplitz,
354 // and thus can not be inverted through the standard Levinson-Durbin method.
355 // There exists a block Levinson-Durbin method, which we may or may not
356 // want to use later. (Eigen's solvers are fast enough that for big matrices,
357 // the convolution operation and not the matrix solving is the bottleneck.)
359 // One thing we definitely want to use, though, is the symmetry properties.
360 // Since we know that g(i, j) = g(|i|, |j|), we can reduce the amount of
361 // unknowns to about 1/4th of the total size. The method is quite simple,
362 // as can be seen from the following toy equation system:
364 // A x0 + B x1 + C x2 = y0
365 // D x0 + E x1 + F x2 = y1
366 // G x0 + H x1 + I x2 = y2
368 // If we now know that e.g. x0=x1 and y0=y1, we can rewrite this to
370 // (A+B+D+E) x0 + (C+F) x2 = 2 y0
371 // (G+H) x0 + I x2 = y2
373 // This both increases accuracy and provides us with a very nice speed
375 MatrixXf M(MatrixXf::Zero((R + 1) * (R + 1), (R + 1) * (R + 1)));
376 MatrixXf r_uv_flattened(MatrixXf::Zero((R + 1) * (R + 1), 1));
377 for (int outer_i = 0; outer_i < 2 * R + 1; ++outer_i) {
378 int folded_outer_i = abs(outer_i - R);
379 for (int outer_j = 0; outer_j < 2 * R + 1; ++outer_j) {
380 int folded_outer_j = abs(outer_j - R);
381 int row = folded_outer_i * (R + 1) + folded_outer_j;
382 for (int inner_i = 0; inner_i < 2 * R + 1; ++inner_i) {
383 int folded_inner_i = abs(inner_i - R);
384 for (int inner_j = 0; inner_j < 2 * R + 1; ++inner_j) {
385 int folded_inner_j = abs(inner_j - R);
386 int col = folded_inner_i * (R + 1) + folded_inner_j;
387 M(row, col) += r_vv((inner_i - R) - (outer_i - R) + 2 * R,
388 (inner_j - R) - (outer_j - R) + 2 * R);
391 r_uv_flattened(row) += r_uv(outer_i, outer_j);
395 LLT<MatrixXf> llt(M);
396 MatrixXf g_flattened = llt.solve(r_uv_flattened);
397 assert(g_flattened.rows() == (R + 1) * (R + 1)),
398 assert(g_flattened.cols() == 1);
400 // Normalize and de-flatten the deconvolution matrix.
401 g = MatrixXf(R + 1, R + 1);
403 for (int i = 0; i < g_flattened.rows(); ++i) {
406 if (y == 0 && x == 0) {
407 sum += g_flattened(i);
408 } else if (y == 0 || x == 0) {
409 sum += 2.0f * g_flattened(i);
411 sum += 4.0f * g_flattened(i);
414 for (int i = 0; i < g_flattened.rows(); ++i) {
417 g(y, x) = g_flattened(i) / sum;
420 last_circle_radius = circle_radius;
421 last_gaussian_radius = gaussian_radius;
422 last_correlation = correlation;
426 void DeconvolutionSharpenEffect::set_gl_state(GLuint glsl_program_num, const std::string &prefix, unsigned *sampler_num)
428 Effect::set_gl_state(glsl_program_num, prefix, sampler_num);
432 if (fabs(circle_radius - last_circle_radius) > 1e-3 ||
433 fabs(gaussian_radius - last_gaussian_radius) > 1e-3 ||
434 fabs(correlation - last_correlation) > 1e-3 ||
435 fabs(noise - last_noise) > 1e-3) {
436 update_deconvolution_kernel();
438 // Now encode it as uniforms, and pass it on to the shader.
439 float samples[4 * (R + 1) * (R + 1)];
440 for (int y = 0; y <= R; ++y) {
441 for (int x = 0; x <= R; ++x) {
442 int i = y * (R + 1) + x;
443 samples[i * 4 + 0] = x / float(width);
444 samples[i * 4 + 1] = y / float(height);
445 samples[i * 4 + 2] = g(y, x);
446 samples[i * 4 + 3] = 0.0f;
450 set_uniform_vec4_array(glsl_program_num, prefix, "samples", samples, (R + 1) * (R + 1));