ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/rt/ambcomp.c
Revision: 2.42
Committed: Wed Apr 30 23:44:06 2014 UTC (10 years ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.41: +7 -6 lines
Log Message:
Minor optimization

File Contents

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