ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/rt/ambcomp.c
Revision: 2.46
Committed: Fri May 2 21:58:50 2014 UTC (10 years ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.45: +197 -61 lines
Log Message:
Added baroque book-keeping to avoid calculating pointless Hessians

File Contents

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