ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/rt/ambcomp.c
Revision: 2.45
Committed: Thu May 1 22:34:25 2014 UTC (10 years ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.44: +7 -7 lines
Log Message:
Fixed bug introduced in change to support -as option (-DNEWAMB)

File Contents

# Content
1 #ifndef lint
2 static const char RCSid[] = "$Id: ambcomp.c,v 2.44 2014/05/01 16:06:11 greg Exp $";
3 #endif
4 /*
5 * Routines to compute "ambient" values using Monte Carlo
6 *
7 * Hessian calculations based on "Practical Hessian-Based Error Control
8 * for Irradiance Caching" by Schwarzhaupt, Wann Jensen, & Jarosz
9 * from ACM SIGGRAPH Asia 2012 conference proceedings.
10 *
11 * Declarations of external symbols in ambient.h
12 */
13
14 #include "copyright.h"
15
16 #include "ray.h"
17 #include "ambient.h"
18 #include "random.h"
19
20 #ifdef NEWAMB
21
22 extern void SDsquare2disk(double ds[2], double seedx, double seedy);
23
24 typedef struct {
25 COLOR v; /* hemisphere sample value */
26 FVECT p; /* intersection point */
27 } AMBSAMP; /* sample value */
28
29 typedef struct {
30 RAY *rp; /* originating ray sample */
31 FVECT ux, uy; /* tangent axis unit vectors */
32 int ns; /* number of samples per axis */
33 COLOR acoef; /* division contribution coefficient */
34 AMBSAMP sa[1]; /* sample array (extends struct) */
35 } AMBHEMI; /* ambient sample hemisphere */
36
37 #define ambsam(h,i,j) (h)->sa[(i)*(h)->ns + (j)]
38
39 typedef struct {
40 FVECT r_i, r_i1, e_i, rcp, rI2_eJ2;
41 double I1, I2;
42 } FFTRI; /* vectors and coefficients for Hessian calculation */
43
44
45 static AMBHEMI *
46 inithemi( /* initialize sampling hemisphere */
47 COLOR ac,
48 RAY *r,
49 double wt
50 )
51 {
52 AMBHEMI *hp;
53 double d;
54 int n, i;
55 /* set number of divisions */
56 if (ambacc <= FTINY &&
57 wt > (d = 0.8*intens(ac)*r->rweight/(ambdiv*minweight)))
58 wt = d; /* avoid ray termination */
59 n = sqrt(ambdiv * wt) + 0.5;
60 i = 1 + 5*(ambacc > FTINY); /* minimum number of samples */
61 if (n < i)
62 n = i;
63 /* allocate sampling array */
64 hp = (AMBHEMI *)malloc(sizeof(AMBHEMI) + sizeof(AMBSAMP)*(n*n - 1));
65 if (hp == NULL)
66 return(NULL);
67 hp->rp = r;
68 hp->ns = n;
69 /* assign coefficient */
70 copycolor(hp->acoef, ac);
71 d = 1.0/(n*n);
72 scalecolor(hp->acoef, d);
73 /* make tangent plane axes */
74 hp->uy[0] = 0.5 - frandom();
75 hp->uy[1] = 0.5 - frandom();
76 hp->uy[2] = 0.5 - frandom();
77 for (i = 3; i--; )
78 if ((-0.6 < r->ron[i]) & (r->ron[i] < 0.6))
79 break;
80 if (i < 0)
81 error(CONSISTENCY, "bad ray direction in inithemi");
82 hp->uy[i] = 1.0;
83 VCROSS(hp->ux, hp->uy, r->ron);
84 normalize(hp->ux);
85 VCROSS(hp->uy, r->ron, hp->ux);
86 /* we're ready to sample */
87 return(hp);
88 }
89
90
91 /* Sample ambient division and apply weighting coefficient */
92 static int
93 getambsamp(RAY *arp, AMBHEMI *hp, int i, int j, int n)
94 {
95 int hlist[3], ii;
96 double spt[2], zd;
97 /* ambient coefficient for weight */
98 if (ambacc > FTINY)
99 setcolor(arp->rcoef, AVGREFL, AVGREFL, AVGREFL);
100 else
101 copycolor(arp->rcoef, hp->acoef);
102 if (rayorigin(arp, AMBIENT, hp->rp, arp->rcoef) < 0)
103 return(0);
104 if (ambacc > FTINY) {
105 multcolor(arp->rcoef, hp->acoef);
106 scalecolor(arp->rcoef, 1./AVGREFL);
107 }
108 hlist[0] = hp->rp->rno;
109 hlist[1] = i;
110 hlist[2] = j;
111 multisamp(spt, 2, urand(ilhash(hlist,3)+n));
112 if (!n) { /* avoid border samples for n==0 */
113 if ((spt[0] < 0.1) | (spt[0] > 0.9))
114 spt[0] = 0.1 + 0.8*frandom();
115 if ((spt[1] < 0.1) | (spt[1] > 0.9))
116 spt[1] = 0.1 + 0.8*frandom();
117 }
118 SDsquare2disk(spt, (i+spt[0])/hp->ns, (j+spt[1])/hp->ns);
119 zd = sqrt(1. - spt[0]*spt[0] - spt[1]*spt[1]);
120 for (ii = 3; ii--; )
121 arp->rdir[ii] = spt[0]*hp->ux[ii] +
122 spt[1]*hp->uy[ii] +
123 zd*hp->rp->ron[ii];
124 checknorm(arp->rdir);
125 dimlist[ndims++] = i*hp->ns + j + 90171;
126 rayvalue(arp); /* evaluate ray */
127 ndims--; /* apply coefficient */
128 multcolor(arp->rcol, arp->rcoef);
129 return(1);
130 }
131
132
133 static AMBSAMP *
134 ambsample( /* initial ambient division sample */
135 AMBHEMI *hp,
136 int i,
137 int j
138 )
139 {
140 AMBSAMP *ap = &ambsam(hp,i,j);
141 RAY ar;
142 /* generate hemispherical sample */
143 if (!getambsamp(&ar, hp, i, j, 0))
144 goto badsample;
145 /* limit vertex distance */
146 if (ar.rt > 10.0*thescene.cusize)
147 ar.rt = 10.0*thescene.cusize;
148 else if (ar.rt <= FTINY) /* should never happen! */
149 goto badsample;
150 VSUM(ap->p, ar.rorg, ar.rdir, ar.rt);
151 copycolor(ap->v, ar.rcol);
152 return(ap);
153 badsample:
154 setcolor(ap->v, 0., 0., 0.);
155 VCOPY(ap->p, hp->rp->rop);
156 return(NULL);
157 }
158
159
160 /* Estimate errors based on ambient division differences */
161 static float *
162 getambdiffs(AMBHEMI *hp)
163 {
164 float *earr = (float *)calloc(hp->ns*hp->ns, sizeof(float));
165 float *ep;
166 AMBSAMP *ap;
167 double b, d2;
168 int i, j;
169
170 if (earr == NULL) /* out of memory? */
171 return(NULL);
172 /* compute squared neighbor diffs */
173 for (ap = hp->sa, ep = earr, i = 0; i < hp->ns; i++)
174 for (j = 0; j < hp->ns; j++, ap++, ep++) {
175 b = bright(ap[0].v);
176 if (i) { /* from above */
177 d2 = b - bright(ap[-hp->ns].v);
178 d2 *= d2;
179 ep[0] += d2;
180 ep[-hp->ns] += d2;
181 }
182 if (j) { /* from behind */
183 d2 = b - bright(ap[-1].v);
184 d2 *= d2;
185 ep[0] += d2;
186 ep[-1] += d2;
187 }
188 }
189 /* correct for number of neighbors */
190 earr[0] *= 2.f;
191 earr[hp->ns-1] *= 2.f;
192 earr[(hp->ns-1)*hp->ns] *= 2.f;
193 earr[(hp->ns-1)*hp->ns + hp->ns-1] *= 2.f;
194 for (i = 1; i < hp->ns-1; i++) {
195 earr[i*hp->ns] *= 4./3.;
196 earr[i*hp->ns + hp->ns-1] *= 4./3.;
197 }
198 for (j = 1; j < hp->ns-1; j++) {
199 earr[j] *= 4./3.;
200 earr[(hp->ns-1)*hp->ns + j] *= 4./3.;
201 }
202 return(earr);
203 }
204
205
206 /* Perform super-sampling on hemisphere (introduces bias) */
207 static void
208 ambsupersamp(double acol[3], AMBHEMI *hp, int cnt)
209 {
210 float *earr = getambdiffs(hp);
211 double e2sum = 0;
212 AMBSAMP *ap;
213 RAY ar;
214 COLOR asum;
215 float *ep;
216 int i, j, n;
217
218 if (earr == NULL) /* just skip calc. if no memory */
219 return;
220 /* add up estimated variances */
221 for (ep = earr + hp->ns*hp->ns; ep-- > earr; )
222 e2sum += *ep;
223 ep = earr; /* perform super-sampling */
224 for (ap = hp->sa, i = 0; i < hp->ns; i++)
225 for (j = 0; j < hp->ns; j++, ap++) {
226 int nss = *ep/e2sum*cnt + frandom();
227 setcolor(asum, 0., 0., 0.);
228 for (n = 1; n <= nss; n++) {
229 if (!getambsamp(&ar, hp, i, j, n)) {
230 nss = n-1;
231 break;
232 }
233 addcolor(asum, ar.rcol);
234 }
235 if (nss) { /* update returned ambient value */
236 const double ssf = 1./(nss + 1);
237 for (n = 3; n--; )
238 acol[n] += ssf*colval(asum,n) +
239 (ssf - 1.)*colval(ap->v,n);
240 }
241 e2sum -= *ep++; /* update remainders */
242 cnt -= nss;
243 }
244 free(earr);
245 }
246
247
248 /* Compute vectors and coefficients for Hessian/gradient calcs */
249 static void
250 comp_fftri(FFTRI *ftp, FVECT ap0, FVECT ap1, FVECT rop)
251 {
252 double rdot_cp, dot_e, dot_er, rdot_r, rdot_r1, J2;
253 int i;
254
255 VSUB(ftp->r_i, ap0, rop);
256 VSUB(ftp->r_i1, ap1, rop);
257 VSUB(ftp->e_i, ap1, ap0);
258 VCROSS(ftp->rcp, ftp->r_i, ftp->r_i1);
259 rdot_cp = 1.0/DOT(ftp->rcp,ftp->rcp);
260 dot_e = DOT(ftp->e_i,ftp->e_i);
261 dot_er = DOT(ftp->e_i, ftp->r_i);
262 rdot_r = 1.0/DOT(ftp->r_i,ftp->r_i);
263 rdot_r1 = 1.0/DOT(ftp->r_i1,ftp->r_i1);
264 ftp->I1 = acos( DOT(ftp->r_i, ftp->r_i1) * sqrt(rdot_r*rdot_r1) ) *
265 sqrt( rdot_cp );
266 ftp->I2 = ( DOT(ftp->e_i, ftp->r_i1)*rdot_r1 - dot_er*rdot_r +
267 dot_e*ftp->I1 )*0.5*rdot_cp;
268 J2 = ( 0.5*(rdot_r - rdot_r1) - dot_er*ftp->I2 ) / dot_e;
269 for (i = 3; i--; )
270 ftp->rI2_eJ2[i] = ftp->I2*ftp->r_i[i] + J2*ftp->e_i[i];
271 }
272
273
274 /* Compose 3x3 matrix from two vectors */
275 static void
276 compose_matrix(FVECT mat[3], FVECT va, FVECT vb)
277 {
278 mat[0][0] = 2.0*va[0]*vb[0];
279 mat[1][1] = 2.0*va[1]*vb[1];
280 mat[2][2] = 2.0*va[2]*vb[2];
281 mat[0][1] = mat[1][0] = va[0]*vb[1] + va[1]*vb[0];
282 mat[0][2] = mat[2][0] = va[0]*vb[2] + va[2]*vb[0];
283 mat[1][2] = mat[2][1] = va[1]*vb[2] + va[2]*vb[1];
284 }
285
286
287 /* Compute partial 3x3 Hessian matrix for edge */
288 static void
289 comp_hessian(FVECT hess[3], FFTRI *ftp, FVECT nrm)
290 {
291 FVECT ncp;
292 FVECT m1[3], m2[3], m3[3], m4[3];
293 double d1, d2, d3, d4;
294 double I3, J3, K3;
295 int i, j;
296 /* compute intermediate coefficients */
297 d1 = 1.0/DOT(ftp->r_i,ftp->r_i);
298 d2 = 1.0/DOT(ftp->r_i1,ftp->r_i1);
299 d3 = 1.0/DOT(ftp->e_i,ftp->e_i);
300 d4 = DOT(ftp->e_i, ftp->r_i);
301 I3 = ( DOT(ftp->e_i, ftp->r_i1)*d2*d2 - d4*d1*d1 + 3.0/d3*ftp->I2 )
302 / ( 4.0*DOT(ftp->rcp,ftp->rcp) );
303 J3 = 0.25*d3*(d1*d1 - d2*d2) - d4*d3*I3;
304 K3 = d3*(ftp->I2 - I3/d1 - 2.0*d4*J3);
305 /* intermediate matrices */
306 VCROSS(ncp, nrm, ftp->e_i);
307 compose_matrix(m1, ncp, ftp->rI2_eJ2);
308 compose_matrix(m2, ftp->r_i, ftp->r_i);
309 compose_matrix(m3, ftp->e_i, ftp->e_i);
310 compose_matrix(m4, ftp->r_i, ftp->e_i);
311 d1 = DOT(nrm, ftp->rcp);
312 d2 = -d1*ftp->I2;
313 d1 *= 2.0;
314 for (i = 3; i--; ) /* final matrix sum */
315 for (j = 3; j--; ) {
316 hess[i][j] = m1[i][j] + d1*( I3*m2[i][j] + K3*m3[i][j] +
317 2.0*J3*m4[i][j] );
318 hess[i][j] += d2*(i==j);
319 hess[i][j] *= 1.0/PI;
320 }
321 }
322
323
324 /* Reverse hessian calculation result for edge in other direction */
325 static void
326 rev_hessian(FVECT hess[3])
327 {
328 int i;
329
330 for (i = 3; i--; ) {
331 hess[i][0] = -hess[i][0];
332 hess[i][1] = -hess[i][1];
333 hess[i][2] = -hess[i][2];
334 }
335 }
336
337
338 /* Add to radiometric Hessian from the given triangle */
339 static void
340 add2hessian(FVECT hess[3], FVECT ehess1[3],
341 FVECT ehess2[3], FVECT ehess3[3], COLORV v)
342 {
343 int i, j;
344
345 for (i = 3; i--; )
346 for (j = 3; j--; )
347 hess[i][j] += v*( ehess1[i][j] + ehess2[i][j] + ehess3[i][j] );
348 }
349
350
351 /* Compute partial displacement form factor gradient for edge */
352 static void
353 comp_gradient(FVECT grad, FFTRI *ftp, FVECT nrm)
354 {
355 FVECT ncp;
356 double f1;
357 int i;
358
359 f1 = 2.0*DOT(nrm, ftp->rcp);
360 VCROSS(ncp, nrm, ftp->e_i);
361 for (i = 3; i--; )
362 grad[i] = (-0.5/PI)*( ftp->I1*ncp[i] + f1*ftp->rI2_eJ2[i] );
363 }
364
365
366 /* Reverse gradient calculation result for edge in other direction */
367 static void
368 rev_gradient(FVECT grad)
369 {
370 grad[0] = -grad[0];
371 grad[1] = -grad[1];
372 grad[2] = -grad[2];
373 }
374
375
376 /* Add to displacement gradient from the given triangle */
377 static void
378 add2gradient(FVECT grad, FVECT egrad1, FVECT egrad2, FVECT egrad3, COLORV v)
379 {
380 int i;
381
382 for (i = 3; i--; )
383 grad[i] += v*( egrad1[i] + egrad2[i] + egrad3[i] );
384 }
385
386
387 /* Return brightness of furthest ambient sample */
388 static COLORV
389 back_ambval(AMBSAMP *ap1, AMBSAMP *ap2, AMBSAMP *ap3, FVECT orig)
390 {
391 COLORV vback;
392 FVECT vec;
393 double d2, d2best;
394
395 VSUB(vec, ap1->p, orig);
396 d2best = DOT(vec,vec);
397 vback = colval(ap1->v,CIEY);
398 VSUB(vec, ap2->p, orig);
399 d2 = DOT(vec,vec);
400 if (d2 > d2best) {
401 d2best = d2;
402 vback = colval(ap2->v,CIEY);
403 }
404 VSUB(vec, ap3->p, orig);
405 d2 = DOT(vec,vec);
406 if (d2 > d2best)
407 return(colval(ap3->v,CIEY));
408 return(vback);
409 }
410
411
412 /* Compute anisotropic radii and eigenvector directions */
413 static int
414 eigenvectors(FVECT uv[2], float ra[2], FVECT hessian[3])
415 {
416 double hess2[2][2];
417 FVECT a, b;
418 double evalue[2], slope1, xmag1;
419 int i;
420 /* project Hessian to sample plane */
421 for (i = 3; i--; ) {
422 a[i] = DOT(hessian[i], uv[0]);
423 b[i] = DOT(hessian[i], uv[1]);
424 }
425 hess2[0][0] = DOT(uv[0], a);
426 hess2[0][1] = DOT(uv[0], b);
427 hess2[1][0] = DOT(uv[1], a);
428 hess2[1][1] = DOT(uv[1], b);
429 /* compute eigenvalue(s) */
430 i = quadratic(evalue, 1.0, -hess2[0][0]-hess2[1][1],
431 hess2[0][0]*hess2[1][1]-hess2[0][1]*hess2[1][0]);
432 if (i == 1) /* double-root (circle) */
433 evalue[1] = evalue[0];
434 if (!i || ((evalue[0] = fabs(evalue[0])) <= FTINY*FTINY) |
435 ((evalue[1] = fabs(evalue[1])) <= FTINY*FTINY) )
436 error(INTERNAL, "bad eigenvalue calculation");
437
438 if (evalue[0] > evalue[1]) {
439 ra[0] = sqrt(sqrt(4.0/evalue[0]));
440 ra[1] = sqrt(sqrt(4.0/evalue[1]));
441 slope1 = evalue[1];
442 } else {
443 ra[0] = sqrt(sqrt(4.0/evalue[1]));
444 ra[1] = sqrt(sqrt(4.0/evalue[0]));
445 slope1 = evalue[0];
446 }
447 /* compute unit eigenvectors */
448 if (fabs(hess2[0][1]) <= FTINY)
449 return; /* uv OK as is */
450 slope1 = (slope1 - hess2[0][0]) / hess2[0][1];
451 xmag1 = sqrt(1.0/(1.0 + slope1*slope1));
452 for (i = 3; i--; ) {
453 b[i] = xmag1*uv[0][i] + slope1*xmag1*uv[1][i];
454 a[i] = slope1*xmag1*uv[0][i] - xmag1*uv[1][i];
455 }
456 VCOPY(uv[0], a);
457 VCOPY(uv[1], b);
458 }
459
460
461 static void
462 ambHessian( /* anisotropic radii & pos. gradient */
463 AMBHEMI *hp,
464 FVECT uv[2], /* returned */
465 float ra[2], /* returned (optional) */
466 float pg[2] /* returned (optional) */
467 )
468 {
469 static char memerrmsg[] = "out of memory in ambHessian()";
470 FVECT (*hessrow)[3] = NULL;
471 FVECT *gradrow = NULL;
472 FVECT hessian[3];
473 FVECT gradient;
474 FFTRI fftr;
475 int i, j;
476 /* be sure to assign unit vectors */
477 VCOPY(uv[0], hp->ux);
478 VCOPY(uv[1], hp->uy);
479 /* clock-wise vertex traversal from sample POV */
480 if (ra != NULL) { /* initialize Hessian row buffer */
481 hessrow = (FVECT (*)[3])malloc(sizeof(FVECT)*3*(hp->ns-1));
482 if (hessrow == NULL)
483 error(SYSTEM, memerrmsg);
484 memset(hessian, 0, sizeof(hessian));
485 } else if (pg == NULL) /* bogus call? */
486 return;
487 if (pg != NULL) { /* initialize form factor row buffer */
488 gradrow = (FVECT *)malloc(sizeof(FVECT)*(hp->ns-1));
489 if (gradrow == NULL)
490 error(SYSTEM, memerrmsg);
491 memset(gradient, 0, sizeof(gradient));
492 }
493 /* compute first row of edges */
494 for (j = 0; j < hp->ns-1; j++) {
495 comp_fftri(&fftr, ambsam(hp,0,j).p,
496 ambsam(hp,0,j+1).p, hp->rp->rop);
497 if (hessrow != NULL)
498 comp_hessian(hessrow[j], &fftr, hp->rp->ron);
499 if (gradrow != NULL)
500 comp_gradient(gradrow[j], &fftr, hp->rp->ron);
501 }
502 /* sum each row of triangles */
503 for (i = 0; i < hp->ns-1; i++) {
504 FVECT hesscol[3]; /* compute first vertical edge */
505 FVECT gradcol;
506 comp_fftri(&fftr, ambsam(hp,i,0).p,
507 ambsam(hp,i+1,0).p, hp->rp->rop);
508 if (hessrow != NULL)
509 comp_hessian(hesscol, &fftr, hp->rp->ron);
510 if (gradrow != NULL)
511 comp_gradient(gradcol, &fftr, hp->rp->ron);
512 for (j = 0; j < hp->ns-1; j++) {
513 FVECT hessdia[3]; /* compute triangle contributions */
514 FVECT graddia;
515 COLORV backg;
516 backg = back_ambval(&ambsam(hp,i,j), &ambsam(hp,i,j+1),
517 &ambsam(hp,i+1,j), hp->rp->rop);
518 /* diagonal (inner) edge */
519 comp_fftri(&fftr, ambsam(hp,i,j+1).p,
520 ambsam(hp,i+1,j).p, hp->rp->rop);
521 if (hessrow != NULL) {
522 comp_hessian(hessdia, &fftr, hp->rp->ron);
523 rev_hessian(hesscol);
524 add2hessian(hessian, hessrow[j], hessdia, hesscol, backg);
525 }
526 if (gradrow != NULL) {
527 comp_gradient(graddia, &fftr, hp->rp->ron);
528 rev_gradient(gradcol);
529 add2gradient(gradient, gradrow[j], graddia, gradcol, backg);
530 }
531 /* initialize edge in next row */
532 comp_fftri(&fftr, ambsam(hp,i+1,j+1).p,
533 ambsam(hp,i+1,j).p, hp->rp->rop);
534 if (hessrow != NULL)
535 comp_hessian(hessrow[j], &fftr, hp->rp->ron);
536 if (gradrow != NULL)
537 comp_gradient(gradrow[j], &fftr, hp->rp->ron);
538 /* new column edge & paired triangle */
539 backg = back_ambval(&ambsam(hp,i,j+1), &ambsam(hp,i+1,j+1),
540 &ambsam(hp,i+1,j), hp->rp->rop);
541 comp_fftri(&fftr, ambsam(hp,i,j+1).p, ambsam(hp,i+1,j+1).p,
542 hp->rp->rop);
543 if (hessrow != NULL) {
544 comp_hessian(hesscol, &fftr, hp->rp->ron);
545 rev_hessian(hessdia);
546 add2hessian(hessian, hessrow[j], hessdia, hesscol, backg);
547 if (i < hp->ns-2)
548 rev_hessian(hessrow[j]);
549 }
550 if (gradrow != NULL) {
551 comp_gradient(gradcol, &fftr, hp->rp->ron);
552 rev_gradient(graddia);
553 add2gradient(gradient, gradrow[j], graddia, gradcol, backg);
554 if (i < hp->ns-2)
555 rev_gradient(gradrow[j]);
556 }
557 }
558 }
559 /* release row buffers */
560 if (hessrow != NULL) free(hessrow);
561 if (gradrow != NULL) free(gradrow);
562
563 if (ra != NULL) /* extract eigenvectors & radii */
564 eigenvectors(uv, ra, hessian);
565 if (pg != NULL) { /* tangential position gradient */
566 pg[0] = DOT(gradient, uv[0]);
567 pg[1] = DOT(gradient, uv[1]);
568 }
569 }
570
571
572 /* Compute direction gradient from a hemispherical sampling */
573 static void
574 ambdirgrad(AMBHEMI *hp, FVECT uv[2], float dg[2])
575 {
576 AMBSAMP *ap;
577 double dgsum[2];
578 int n;
579 FVECT vd;
580 double gfact;
581
582 dgsum[0] = dgsum[1] = 0.0; /* sum values times -tan(theta) */
583 for (ap = hp->sa, n = hp->ns*hp->ns; n--; ap++) {
584 /* use vector for azimuth + 90deg */
585 VSUB(vd, ap->p, hp->rp->rop);
586 /* brightness over cosine factor */
587 gfact = colval(ap->v,CIEY) / DOT(hp->rp->ron, vd);
588 /* sine = proj_radius/vd_length */
589 dgsum[0] -= DOT(uv[1], vd) * gfact;
590 dgsum[1] += DOT(uv[0], vd) * gfact;
591 }
592 dg[0] = dgsum[0] / (hp->ns*hp->ns);
593 dg[1] = dgsum[1] / (hp->ns*hp->ns);
594 }
595
596
597 int
598 doambient( /* compute ambient component */
599 COLOR rcol, /* input/output color */
600 RAY *r,
601 double wt,
602 FVECT uv[2], /* returned (optional) */
603 float ra[2], /* returned (optional) */
604 float pg[2], /* returned (optional) */
605 float dg[2] /* returned (optional) */
606 )
607 {
608 AMBHEMI *hp = inithemi(rcol, r, wt);
609 int cnt;
610 FVECT my_uv[2];
611 double d, K, acol[3];
612 AMBSAMP *ap;
613 int i, j;
614 /* check/initialize */
615 if (hp == NULL)
616 return(0);
617 if (uv != NULL)
618 memset(uv, 0, sizeof(FVECT)*2);
619 if (ra != NULL)
620 ra[0] = ra[1] = 0.0;
621 if (pg != NULL)
622 pg[0] = pg[1] = 0.0;
623 if (dg != NULL)
624 dg[0] = dg[1] = 0.0;
625 /* sample the hemisphere */
626 acol[0] = acol[1] = acol[2] = 0.0;
627 cnt = 0;
628 for (i = hp->ns; i--; )
629 for (j = hp->ns; j--; )
630 if ((ap = ambsample(hp, i, j)) != NULL) {
631 addcolor(acol, ap->v);
632 ++cnt;
633 }
634 if (!cnt) {
635 setcolor(rcol, 0.0, 0.0, 0.0);
636 free(hp);
637 return(0); /* no valid samples */
638 }
639 if (cnt < hp->ns*hp->ns) { /* incomplete sampling? */
640 copycolor(rcol, acol);
641 free(hp);
642 return(-1); /* return value w/o Hessian */
643 }
644 cnt = ambssamp*wt + 0.5; /* perform super-sampling? */
645 if (cnt > 0)
646 ambsupersamp(acol, hp, cnt);
647 copycolor(rcol, acol); /* final indirect irradiance/PI */
648 if ((ra == NULL) & (pg == NULL) & (dg == NULL)) {
649 free(hp);
650 return(-1); /* no radius or gradient calc. */
651 }
652 if ((d = bright(acol)) > FTINY) { /* normalize Y values */
653 d = 0.99*(hp->ns*hp->ns)/d;
654 K = 0.01;
655 } else { /* or fall back on geometric Hessian */
656 K = 1.0;
657 pg = NULL;
658 dg = NULL;
659 }
660 ap = hp->sa; /* relative Y channel from here on... */
661 for (i = hp->ns*hp->ns; i--; ap++)
662 colval(ap->v,CIEY) = bright(ap->v)*d + K;
663
664 if (uv == NULL) /* make sure we have axis pointers */
665 uv = my_uv;
666 /* compute radii & pos. gradient */
667 ambHessian(hp, uv, ra, pg);
668
669 if (dg != NULL) /* compute direction gradient */
670 ambdirgrad(hp, uv, dg);
671
672 if (ra != NULL) { /* scale/clamp radii */
673 if (pg != NULL) {
674 if (ra[0]*(d = fabs(pg[0])) > 1.0)
675 ra[0] = 1.0/d;
676 if (ra[1]*(d = fabs(pg[1])) > 1.0)
677 ra[1] = 1.0/d;
678 if (ra[0] > ra[1])
679 ra[0] = ra[1];
680 }
681 if (ra[0] < minarad) {
682 ra[0] = minarad;
683 if (ra[1] < minarad)
684 ra[1] = minarad;
685 }
686 ra[0] *= d = 1.0/sqrt(sqrt(wt));
687 if ((ra[1] *= d) > 2.0*ra[0])
688 ra[1] = 2.0*ra[0];
689 if (ra[1] > maxarad) {
690 ra[1] = maxarad;
691 if (ra[0] > maxarad)
692 ra[0] = maxarad;
693 }
694 if (pg != NULL) { /* cap gradient if necessary */
695 d = pg[0]*pg[0]*ra[0]*ra[0] + pg[1]*pg[1]*ra[1]*ra[1];
696 if (d > 1.0) {
697 d = 1.0/sqrt(d);
698 pg[0] *= d;
699 pg[1] *= d;
700 }
701 }
702 }
703 free(hp); /* clean up and return */
704 return(1);
705 }
706
707
708 #else /* ! NEWAMB */
709
710
711 void
712 inithemi( /* initialize sampling hemisphere */
713 AMBHEMI *hp,
714 COLOR ac,
715 RAY *r,
716 double wt
717 )
718 {
719 double d;
720 int i;
721 /* set number of divisions */
722 if (ambacc <= FTINY &&
723 wt > (d = 0.8*intens(ac)*r->rweight/(ambdiv*minweight)))
724 wt = d; /* avoid ray termination */
725 hp->nt = sqrt(ambdiv * wt / PI) + 0.5;
726 i = ambacc > FTINY ? 3 : 1; /* minimum number of samples */
727 if (hp->nt < i)
728 hp->nt = i;
729 hp->np = PI * hp->nt + 0.5;
730 /* set number of super-samples */
731 hp->ns = ambssamp * wt + 0.5;
732 /* assign coefficient */
733 copycolor(hp->acoef, ac);
734 d = 1.0/(hp->nt*hp->np);
735 scalecolor(hp->acoef, d);
736 /* make axes */
737 VCOPY(hp->uz, r->ron);
738 hp->uy[0] = hp->uy[1] = hp->uy[2] = 0.0;
739 for (i = 0; i < 3; i++)
740 if (hp->uz[i] < 0.6 && hp->uz[i] > -0.6)
741 break;
742 if (i >= 3)
743 error(CONSISTENCY, "bad ray direction in inithemi");
744 hp->uy[i] = 1.0;
745 fcross(hp->ux, hp->uy, hp->uz);
746 normalize(hp->ux);
747 fcross(hp->uy, hp->uz, hp->ux);
748 }
749
750
751 int
752 divsample( /* sample a division */
753 AMBSAMP *dp,
754 AMBHEMI *h,
755 RAY *r
756 )
757 {
758 RAY ar;
759 int hlist[3];
760 double spt[2];
761 double xd, yd, zd;
762 double b2;
763 double phi;
764 int i;
765 /* ambient coefficient for weight */
766 if (ambacc > FTINY)
767 setcolor(ar.rcoef, AVGREFL, AVGREFL, AVGREFL);
768 else
769 copycolor(ar.rcoef, h->acoef);
770 if (rayorigin(&ar, AMBIENT, r, ar.rcoef) < 0)
771 return(-1);
772 if (ambacc > FTINY) {
773 multcolor(ar.rcoef, h->acoef);
774 scalecolor(ar.rcoef, 1./AVGREFL);
775 }
776 hlist[0] = r->rno;
777 hlist[1] = dp->t;
778 hlist[2] = dp->p;
779 multisamp(spt, 2, urand(ilhash(hlist,3)+dp->n));
780 zd = sqrt((dp->t + spt[0])/h->nt);
781 phi = 2.0*PI * (dp->p + spt[1])/h->np;
782 xd = tcos(phi) * zd;
783 yd = tsin(phi) * zd;
784 zd = sqrt(1.0 - zd*zd);
785 for (i = 0; i < 3; i++)
786 ar.rdir[i] = xd*h->ux[i] +
787 yd*h->uy[i] +
788 zd*h->uz[i];
789 checknorm(ar.rdir);
790 dimlist[ndims++] = dp->t*h->np + dp->p + 90171;
791 rayvalue(&ar);
792 ndims--;
793 multcolor(ar.rcol, ar.rcoef); /* apply coefficient */
794 addcolor(dp->v, ar.rcol);
795 /* use rt to improve gradient calc */
796 if (ar.rt > FTINY && ar.rt < FHUGE)
797 dp->r += 1.0/ar.rt;
798 /* (re)initialize error */
799 if (dp->n++) {
800 b2 = bright(dp->v)/dp->n - bright(ar.rcol);
801 b2 = b2*b2 + dp->k*((dp->n-1)*(dp->n-1));
802 dp->k = b2/(dp->n*dp->n);
803 } else
804 dp->k = 0.0;
805 return(0);
806 }
807
808
809 static int
810 ambcmp( /* decreasing order */
811 const void *p1,
812 const void *p2
813 )
814 {
815 const AMBSAMP *d1 = (const AMBSAMP *)p1;
816 const AMBSAMP *d2 = (const AMBSAMP *)p2;
817
818 if (d1->k < d2->k)
819 return(1);
820 if (d1->k > d2->k)
821 return(-1);
822 return(0);
823 }
824
825
826 static int
827 ambnorm( /* standard order */
828 const void *p1,
829 const void *p2
830 )
831 {
832 const AMBSAMP *d1 = (const AMBSAMP *)p1;
833 const AMBSAMP *d2 = (const AMBSAMP *)p2;
834 int c;
835
836 if ( (c = d1->t - d2->t) )
837 return(c);
838 return(d1->p - d2->p);
839 }
840
841
842 double
843 doambient( /* compute ambient component */
844 COLOR rcol,
845 RAY *r,
846 double wt,
847 FVECT pg,
848 FVECT dg
849 )
850 {
851 double b, d=0;
852 AMBHEMI hemi;
853 AMBSAMP *div;
854 AMBSAMP dnew;
855 double acol[3];
856 AMBSAMP *dp;
857 double arad;
858 int divcnt;
859 int i, j;
860 /* initialize hemisphere */
861 inithemi(&hemi, rcol, r, wt);
862 divcnt = hemi.nt * hemi.np;
863 /* initialize */
864 if (pg != NULL)
865 pg[0] = pg[1] = pg[2] = 0.0;
866 if (dg != NULL)
867 dg[0] = dg[1] = dg[2] = 0.0;
868 setcolor(rcol, 0.0, 0.0, 0.0);
869 if (divcnt == 0)
870 return(0.0);
871 /* allocate super-samples */
872 if (hemi.ns > 0 || pg != NULL || dg != NULL) {
873 div = (AMBSAMP *)malloc(divcnt*sizeof(AMBSAMP));
874 if (div == NULL)
875 error(SYSTEM, "out of memory in doambient");
876 } else
877 div = NULL;
878 /* sample the divisions */
879 arad = 0.0;
880 acol[0] = acol[1] = acol[2] = 0.0;
881 if ((dp = div) == NULL)
882 dp = &dnew;
883 divcnt = 0;
884 for (i = 0; i < hemi.nt; i++)
885 for (j = 0; j < hemi.np; j++) {
886 dp->t = i; dp->p = j;
887 setcolor(dp->v, 0.0, 0.0, 0.0);
888 dp->r = 0.0;
889 dp->n = 0;
890 if (divsample(dp, &hemi, r) < 0) {
891 if (div != NULL)
892 dp++;
893 continue;
894 }
895 arad += dp->r;
896 divcnt++;
897 if (div != NULL)
898 dp++;
899 else
900 addcolor(acol, dp->v);
901 }
902 if (!divcnt) {
903 if (div != NULL)
904 free((void *)div);
905 return(0.0); /* no samples taken */
906 }
907 if (divcnt < hemi.nt*hemi.np) {
908 pg = dg = NULL; /* incomplete sampling */
909 hemi.ns = 0;
910 } else if (arad > FTINY && divcnt/arad < minarad) {
911 hemi.ns = 0; /* close enough */
912 } else if (hemi.ns > 0) { /* else perform super-sampling? */
913 comperrs(div, &hemi); /* compute errors */
914 qsort(div, divcnt, sizeof(AMBSAMP), ambcmp); /* sort divs */
915 /* super-sample */
916 for (i = hemi.ns; i > 0; i--) {
917 dnew = *div;
918 if (divsample(&dnew, &hemi, r) < 0) {
919 dp++;
920 continue;
921 }
922 dp = div; /* reinsert */
923 j = divcnt < i ? divcnt : i;
924 while (--j > 0 && dnew.k < dp[1].k) {
925 *dp = *(dp+1);
926 dp++;
927 }
928 *dp = dnew;
929 }
930 if (pg != NULL || dg != NULL) /* restore order */
931 qsort(div, divcnt, sizeof(AMBSAMP), ambnorm);
932 }
933 /* compute returned values */
934 if (div != NULL) {
935 arad = 0.0; /* note: divcnt may be < nt*np */
936 for (i = hemi.nt*hemi.np, dp = div; i-- > 0; dp++) {
937 arad += dp->r;
938 if (dp->n > 1) {
939 b = 1.0/dp->n;
940 scalecolor(dp->v, b);
941 dp->r *= b;
942 dp->n = 1;
943 }
944 addcolor(acol, dp->v);
945 }
946 b = bright(acol);
947 if (b > FTINY) {
948 b = 1.0/b; /* compute & normalize gradient(s) */
949 if (pg != NULL) {
950 posgradient(pg, div, &hemi);
951 for (i = 0; i < 3; i++)
952 pg[i] *= b;
953 }
954 if (dg != NULL) {
955 dirgradient(dg, div, &hemi);
956 for (i = 0; i < 3; i++)
957 dg[i] *= b;
958 }
959 }
960 free((void *)div);
961 }
962 copycolor(rcol, acol);
963 if (arad <= FTINY)
964 arad = maxarad;
965 else
966 arad = (divcnt+hemi.ns)/arad;
967 if (pg != NULL) { /* reduce radius if gradient large */
968 d = DOT(pg,pg);
969 if (d*arad*arad > 1.0)
970 arad = 1.0/sqrt(d);
971 }
972 if (arad < minarad) {
973 arad = minarad;
974 if (pg != NULL && d*arad*arad > 1.0) { /* cap gradient */
975 d = 1.0/arad/sqrt(d);
976 for (i = 0; i < 3; i++)
977 pg[i] *= d;
978 }
979 }
980 if ((arad /= sqrt(wt)) > maxarad)
981 arad = maxarad;
982 return(arad);
983 }
984
985
986 void
987 comperrs( /* compute initial error estimates */
988 AMBSAMP *da, /* assumes standard ordering */
989 AMBHEMI *hp
990 )
991 {
992 double b, b2;
993 int i, j;
994 AMBSAMP *dp;
995 /* sum differences from neighbors */
996 dp = da;
997 for (i = 0; i < hp->nt; i++)
998 for (j = 0; j < hp->np; j++) {
999 #ifdef DEBUG
1000 if (dp->t != i || dp->p != j)
1001 error(CONSISTENCY,
1002 "division order in comperrs");
1003 #endif
1004 b = bright(dp[0].v);
1005 if (i > 0) { /* from above */
1006 b2 = bright(dp[-hp->np].v) - b;
1007 b2 *= b2 * 0.25;
1008 dp[0].k += b2;
1009 dp[-hp->np].k += b2;
1010 }
1011 if (j > 0) { /* from behind */
1012 b2 = bright(dp[-1].v) - b;
1013 b2 *= b2 * 0.25;
1014 dp[0].k += b2;
1015 dp[-1].k += b2;
1016 } else { /* around */
1017 b2 = bright(dp[hp->np-1].v) - b;
1018 b2 *= b2 * 0.25;
1019 dp[0].k += b2;
1020 dp[hp->np-1].k += b2;
1021 }
1022 dp++;
1023 }
1024 /* divide by number of neighbors */
1025 dp = da;
1026 for (j = 0; j < hp->np; j++) /* top row */
1027 (dp++)->k *= 1.0/3.0;
1028 if (hp->nt < 2)
1029 return;
1030 for (i = 1; i < hp->nt-1; i++) /* central region */
1031 for (j = 0; j < hp->np; j++)
1032 (dp++)->k *= 0.25;
1033 for (j = 0; j < hp->np; j++) /* bottom row */
1034 (dp++)->k *= 1.0/3.0;
1035 }
1036
1037
1038 void
1039 posgradient( /* compute position gradient */
1040 FVECT gv,
1041 AMBSAMP *da, /* assumes standard ordering */
1042 AMBHEMI *hp
1043 )
1044 {
1045 int i, j;
1046 double nextsine, lastsine, b, d;
1047 double mag0, mag1;
1048 double phi, cosp, sinp, xd, yd;
1049 AMBSAMP *dp;
1050
1051 xd = yd = 0.0;
1052 for (j = 0; j < hp->np; j++) {
1053 dp = da + j;
1054 mag0 = mag1 = 0.0;
1055 lastsine = 0.0;
1056 for (i = 0; i < hp->nt; i++) {
1057 #ifdef DEBUG
1058 if (dp->t != i || dp->p != j)
1059 error(CONSISTENCY,
1060 "division order in posgradient");
1061 #endif
1062 b = bright(dp->v);
1063 if (i > 0) {
1064 d = dp[-hp->np].r;
1065 if (dp[0].r > d) d = dp[0].r;
1066 /* sin(t)*cos(t)^2 */
1067 d *= lastsine * (1.0 - (double)i/hp->nt);
1068 mag0 += d*(b - bright(dp[-hp->np].v));
1069 }
1070 nextsine = sqrt((double)(i+1)/hp->nt);
1071 if (j > 0) {
1072 d = dp[-1].r;
1073 if (dp[0].r > d) d = dp[0].r;
1074 mag1 += d * (nextsine - lastsine) *
1075 (b - bright(dp[-1].v));
1076 } else {
1077 d = dp[hp->np-1].r;
1078 if (dp[0].r > d) d = dp[0].r;
1079 mag1 += d * (nextsine - lastsine) *
1080 (b - bright(dp[hp->np-1].v));
1081 }
1082 dp += hp->np;
1083 lastsine = nextsine;
1084 }
1085 mag0 *= 2.0*PI / hp->np;
1086 phi = 2.0*PI * (double)j/hp->np;
1087 cosp = tcos(phi); sinp = tsin(phi);
1088 xd += mag0*cosp - mag1*sinp;
1089 yd += mag0*sinp + mag1*cosp;
1090 }
1091 for (i = 0; i < 3; i++)
1092 gv[i] = (xd*hp->ux[i] + yd*hp->uy[i])*(hp->nt*hp->np)/PI;
1093 }
1094
1095
1096 void
1097 dirgradient( /* compute direction gradient */
1098 FVECT gv,
1099 AMBSAMP *da, /* assumes standard ordering */
1100 AMBHEMI *hp
1101 )
1102 {
1103 int i, j;
1104 double mag;
1105 double phi, xd, yd;
1106 AMBSAMP *dp;
1107
1108 xd = yd = 0.0;
1109 for (j = 0; j < hp->np; j++) {
1110 dp = da + j;
1111 mag = 0.0;
1112 for (i = 0; i < hp->nt; i++) {
1113 #ifdef DEBUG
1114 if (dp->t != i || dp->p != j)
1115 error(CONSISTENCY,
1116 "division order in dirgradient");
1117 #endif
1118 /* tan(t) */
1119 mag += bright(dp->v)/sqrt(hp->nt/(i+.5) - 1.0);
1120 dp += hp->np;
1121 }
1122 phi = 2.0*PI * (j+.5)/hp->np + PI/2.0;
1123 xd += mag * tcos(phi);
1124 yd += mag * tsin(phi);
1125 }
1126 for (i = 0; i < 3; i++)
1127 gv[i] = xd*hp->ux[i] + yd*hp->uy[i];
1128 }
1129
1130 #endif /* ! NEWAMB */