ViewVC Help
View File | Revision Log | Show Annotations | Download File | Root Listing
root/radiance/ray/src/util/rfluxmtx.c
Revision: 2.8
Committed: Fri Jul 25 15:46:10 2014 UTC (9 years, 8 months ago) by greg
Content type: text/plain
Branch: MAIN
Changes since 2.7: +20 -5 lines
Log Message:
Added -I option support in pass-through and fixed source orientation

File Contents

# User Rev Content
1 greg 2.1 #ifndef lint
2 greg 2.8 static const char RCSid[] = "$Id: rfluxmtx.c,v 2.7 2014/07/23 22:39:40 greg Exp $";
3 greg 2.1 #endif
4     /*
5     * Calculate flux transfer matrix or matrices using rcontrib
6     */
7    
8     #include "copyright.h"
9    
10     #include <ctype.h>
11     #include <stdlib.h>
12     #include "rtio.h"
13     #include "rtmath.h"
14     #include "bsdf.h"
15     #include "bsdf_m.h"
16     #include "random.h"
17     #include "triangulate.h"
18     #include "platform.h"
19    
20     #ifdef getc_unlocked /* avoid horrendous overhead of flockfile */
21     #undef getc
22     #define getc getc_unlocked
23     #endif
24    
25     #ifdef _WIN32
26 greg 2.2 #define SPECIALS " \t\"$*?"
27 greg 2.1 #define QUOTCHAR '"'
28     #else
29 greg 2.2 #define SPECIALS " \t\n'\"()${}*?[];|&"
30 greg 2.1 #define QUOTCHAR '\''
31     #define ALTQUOT '"'
32     #endif
33    
34     #define MAXRCARG 512
35    
36     char *progname; /* global argv[0] */
37    
38     int verbose = 0; /* verbose mode? */
39    
40     char *rcarg[MAXRCARG+1] = {"rcontrib", "-fo+"};
41     int nrcargs = 2;
42    
43     const char overflowerr[] = "%s: too many arguments!\n";
44    
45     #define CHECKARGC(n) if (nrcargs >= MAXRCARG-(n)) \
46     { fprintf(stderr, overflowerr, progname); exit(1); }
47    
48 greg 2.2 int sampcnt = 0; /* sample count (0==unset) */
49 greg 2.1
50 greg 2.2 char *reinhfn = "reinhartb.cal";
51     char *shirchiufn = "disk2square.cal";
52 greg 2.1 char *kfullfn = "klems_full.cal";
53     char *khalffn = "klems_half.cal";
54     char *kquarterfn = "klems_quarter.cal";
55    
56 greg 2.5 /* string indicating parameters */
57     const char PARAMSTART[] = "@rfluxmtx";
58 greg 2.1
59     /* surface type IDs */
60     #define ST_NONE 0
61     #define ST_POLY 1
62     #define ST_RING 2
63     #define ST_SOURCE 3
64    
65     typedef struct surf_s {
66     struct surf_s *next; /* next surface in list */
67     void *priv; /* private data (malloc'ed) */
68     char sname[32]; /* surface name */
69     FVECT snrm; /* surface normal */
70 greg 2.4 double area; /* surface area / proj. solid angle */
71 greg 2.1 short styp; /* surface type */
72 greg 2.2 short nfargs; /* number of real arguments */
73     double farg[1]; /* real values (extends struct) */
74 greg 2.1 } SURF; /* surface structure */
75    
76     typedef struct {
77     FVECT uva[2]; /* tangent axes */
78     int ntris; /* number of triangles */
79     struct ptri {
80     float afrac; /* fraction of total area */
81     short vndx[3]; /* vertex indices */
82     } tri[1]; /* triangle array (extends struct) */
83     } POLYTRIS; /* triangulated polygon */
84    
85     typedef struct param_s {
86     char hemis[32]; /* hemispherical sampling spec. */
87     int hsiz; /* hemisphere basis size */
88     int nsurfs; /* number of surfaces */
89     SURF *slist; /* list of surfaces */
90     FVECT vup; /* up vector (zero if unset) */
91     FVECT nrm; /* average normal direction */
92     FVECT udir, vdir; /* v-up tangent axes */
93     char *outfn; /* output file name (receiver) */
94     int (*sample_basis)(struct param_s *p, int, FILE *);
95     } PARAMS; /* sender/receiver parameters */
96    
97     PARAMS curparams;
98     char curmod[128];
99 greg 2.6 char newparams[1024];
100 greg 2.1
101     typedef int SURFSAMP(FVECT, SURF *, double);
102    
103     static SURFSAMP ssamp_bad, ssamp_poly, ssamp_ring;
104    
105     SURFSAMP *orig_in_surf[4] = {
106     ssamp_bad, ssamp_poly, ssamp_ring, ssamp_bad
107     };
108    
109     /* Clear parameter set */
110     static void
111     clear_params(PARAMS *p, int reset_only)
112     {
113     while (p->slist != NULL) {
114     SURF *sdel = p->slist;
115     p->slist = sdel->next;
116     if (sdel->priv != NULL)
117     free(sdel->priv);
118     free(sdel);
119     }
120     if (reset_only) {
121     p->nsurfs = 0;
122     memset(p->nrm, 0, sizeof(FVECT));
123     memset(p->vup, 0, sizeof(FVECT));
124     p->outfn = NULL;
125     return;
126     }
127 greg 2.6 memset(p, 0, sizeof(PARAMS));
128 greg 2.1 }
129    
130     /* Get surface type from name */
131     static int
132     surf_type(const char *otype)
133     {
134     if (!strcmp(otype, "polygon"))
135     return(ST_POLY);
136     if (!strcmp(otype, "ring"))
137     return(ST_RING);
138     if (!strcmp(otype, "source"))
139     return(ST_SOURCE);
140     return(ST_NONE);
141     }
142    
143     /* Add arguments to oconv command */
144     static char *
145     oconv_command(int ac, char *av[])
146     {
147     static char oconvbuf[2048] = "!oconv -f";
148     char *cp = oconvbuf + 9;
149    
150     while (ac-- > 0) {
151     if (cp >= oconvbuf+(sizeof(oconvbuf)-32)) {
152     fputs(progname, stderr);
153     fputs(": too many file arguments!\n", stderr);
154     exit(1);
155     }
156     *cp++ = ' ';
157     strcpy(cp, *av++);
158     while (*cp) cp++;
159     }
160     *cp = '\0';
161     return(oconvbuf);
162     }
163    
164     /* Check if any of the characters in str2 are found in str1 */
165     static int
166     matchany(const char *str1, const char *str2)
167     {
168     while (*str1) {
169     const char *cp = str2;
170     while (*cp)
171     if (*cp++ == *str1)
172     return(*str1);
173     ++str1;
174     }
175     return(0);
176     }
177    
178    
179     /* Convert a set of arguments into a command line for pipe() or system() */
180     static char *
181     convert_commandline(char *cmd, const int len, char *av[])
182     {
183     int match;
184     char *cp;
185    
186     for (cp = cmd; *av != NULL; av++) {
187     const int n = strlen(*av);
188     if (cp+n >= cmd+(len-3)) {
189     fputs(progname, stderr);
190     return(NULL);
191     }
192     if ((match = matchany(*av, SPECIALS))) {
193     const int quote =
194     #ifdef ALTQUOT
195     (match == QUOTCHAR) ? ALTQUOT :
196     #endif
197     QUOTCHAR;
198     *cp++ = quote;
199     strcpy(cp, *av);
200     cp += n;
201     *cp++ = quote;
202     } else {
203     strcpy(cp, *av);
204     cp += n;
205     }
206     *cp++ = ' ';
207     }
208     if (cp <= cmd)
209     return(NULL);
210     *--cp = '\0';
211     return(cmd);
212     }
213    
214     /* Open a pipe to/from a command given as an argument list */
215     static FILE *
216     popen_arglist(char *av[], char *mode)
217     {
218     char cmd[10240];
219    
220     if (!convert_commandline(cmd, sizeof(cmd), av)) {
221     fputs(progname, stderr);
222     fputs(": command line too long in popen_arglist()\n", stderr);
223     return(NULL);
224     }
225     if (verbose)
226     fprintf(stderr, "%s: opening pipe %s: %s\n",
227     progname, (*mode=='w') ? "to" : "from", cmd);
228     return(popen(cmd, mode));
229     }
230    
231     #ifdef _WIN32
232     /* Execute system command (Windows version) */
233     static int
234     my_exec(char *av[])
235     {
236     char cmd[10240];
237    
238     if (!convert_commandline(cmd, sizeof(cmd), av)) {
239     fputs(progname, stderr);
240     fputs(": command line too long in my_exec()\n", stderr);
241     return(1);
242     }
243     if (verbose)
244     fprintf(stderr, "%s: running: %s\n", progname, cmd);
245     return(system(cmd));
246     }
247     #else
248     /* Execute system command in our stead (Unix version) */
249     static int
250     my_exec(char *av[])
251     {
252     char *compath;
253    
254     if ((compath = getpath((char *)av[0], getenv("PATH"), X_OK)) == NULL) {
255     fprintf(stderr, "%s: cannot locate %s\n", progname, av[0]);
256     return(1);
257     }
258     if (verbose) {
259     char cmd[4096];
260     if (!convert_commandline(cmd, sizeof(cmd), av))
261     strcpy(cmd, "COMMAND TOO LONG TO SHOW");
262     fprintf(stderr, "%s: running: %s\n", progname, cmd);
263     }
264     execv(compath, av); /* successful call never returns */
265     perror(compath);
266     return(1);
267     }
268     #endif
269    
270     /* Get normalized direction vector from string specification */
271     static int
272     get_direction(FVECT dv, const char *s)
273     {
274     int sign = 1;
275     int axis = 0;
276    
277     memset(dv, 0, sizeof(FVECT));
278     nextchar:
279     switch (*s) {
280     case '+':
281     ++s;
282     goto nextchar;
283     case '-':
284     sign = -sign;
285     ++s;
286     goto nextchar;
287     case 'z':
288     case 'Z':
289     ++axis;
290     case 'y':
291     case 'Y':
292     ++axis;
293     case 'x':
294     case 'X':
295     dv[axis] = sign;
296     return(!s[1] | isspace(s[1]));
297     default:
298     break;
299     }
300     #ifdef SMLFLT
301     if (sscanf(s, "%f,%f,%f", &dv[0], &dv[1], &dv[2]) != 3)
302     #else
303     if (sscanf(s, "%lf,%lf,%lf", &dv[0], &dv[1], &dv[2]) != 3)
304     #endif
305     return(0);
306     dv[0] *= (RREAL)sign;
307     return(normalize(dv) > 0);
308     }
309    
310     /* Parse program parameters (directives) */
311     static int
312 greg 2.6 parse_params(PARAMS *p, char *pargs)
313 greg 2.1 {
314     char *cp = pargs;
315     int nparams = 0;
316     int i;
317    
318     for ( ; ; )
319     switch (*cp++) {
320     case 'h':
321     if (*cp++ != '=')
322     break;
323 greg 2.6 p->hsiz = 0;
324 greg 2.1 i = 0;
325     while (*cp && !isspace(*cp)) {
326     if (isdigit(*cp))
327 greg 2.6 p->hsiz = 10*p->hsiz + *cp - '0';
328     p->hemis[i++] = *cp++;
329 greg 2.1 }
330     if (!i)
331     break;
332 greg 2.6 p->hemis[i] = '\0';
333     p->hsiz += !p->hsiz;
334 greg 2.1 ++nparams;
335     continue;
336     case 'u':
337     if (*cp++ != '=')
338     break;
339 greg 2.6 if (!get_direction(p->vup, cp))
340 greg 2.1 break;
341     ++nparams;
342     continue;
343     case 'o':
344     if (*cp++ != '=')
345     break;
346     i = 0;
347     while (*cp && !isspace(*cp++))
348     i++;
349     if (!i)
350     break;
351     *--cp = '\0';
352 greg 2.6 p->outfn = savqstr(cp-i);
353 greg 2.1 *cp++ = ' ';
354     ++nparams;
355     continue;
356     case ' ':
357     case '\t':
358     case '\r':
359     case '\n':
360     continue;
361     case '\0':
362     return(nparams);
363     default:
364     break;
365     }
366     fprintf(stderr, "%s: bad parameter string '%s'\n", progname, pargs);
367     exit(1);
368     return(-1); /* pro forma return */
369     }
370    
371     /* Add receiver arguments (directives) corresponding to the current modifier */
372     static void
373     finish_receiver(void)
374     {
375 greg 2.4 char sbuf[256];
376     int uniform = 0;
377 greg 2.2 char *calfn = NULL;
378 greg 2.3 char *params = NULL;
379 greg 2.2 char *binv = NULL;
380     char *binf = NULL;
381     char *nbins = NULL;
382 greg 2.1
383     if (!curmod[0]) {
384     fputs(progname, stderr);
385     fputs(": missing receiver surface!\n", stderr);
386     exit(1);
387     }
388     if (curparams.outfn != NULL) { /* add output file spec. */
389     CHECKARGC(2);
390     rcarg[nrcargs++] = "-o";
391     rcarg[nrcargs++] = curparams.outfn;
392     }
393 greg 2.4 /* check arguments */
394 greg 2.1 if (!curparams.hemis[0]) {
395     fputs(progname, stderr);
396     fputs(": missing hemisphere sampling type!\n", stderr);
397     exit(1);
398     }
399     if (normalize(curparams.nrm) == 0) {
400     fputs(progname, stderr);
401     fputs(": undefined normal for hemisphere sampling\n", stderr);
402     exit(1);
403     }
404     if (normalize(curparams.vup) == 0)
405     if (fabs(curparams.nrm[2]) < .7)
406     curparams.vup[2] = 1;
407     else
408     curparams.vup[1] = 1;
409 greg 2.4 /* determine sample type/bin */
410 greg 2.1 if (tolower(curparams.hemis[0]) == 'u' | curparams.hemis[0] == '1') {
411 greg 2.3 binv = "0"; /* uniform sampling -- one bin */
412 greg 2.4 uniform = 1;
413 greg 2.1 } else if (tolower(curparams.hemis[0]) == 's' &&
414     tolower(curparams.hemis[1]) == 'c') {
415 greg 2.2 /* assign parameters */
416 greg 2.1 if (curparams.hsiz <= 1) {
417     fputs(progname, stderr);
418     fputs(": missing size for Shirley-Chiu sampling!\n", stderr);
419     exit(1);
420     }
421 greg 2.3 calfn = shirchiufn; shirchiufn = NULL;
422 greg 2.2 sprintf(sbuf, "SCdim=%d,Nx=%g,Ny=%g,Nz=%g,Ux=%g,Uy=%g,Uz=%g",
423     curparams.hsiz,
424     curparams.nrm[0], curparams.nrm[1], curparams.nrm[2],
425     curparams.vup[0], curparams.vup[1], curparams.vup[2]);
426 greg 2.3 params = savqstr(sbuf);
427 greg 2.2 binv = "scbin";
428     nbins = "SCdim*SCdim";
429 greg 2.1 } else if ((tolower(curparams.hemis[0]) == 'r') |
430     (tolower(curparams.hemis[0]) == 't')) {
431 greg 2.2 calfn = reinhfn; reinhfn = NULL;
432 greg 2.3 sprintf(sbuf, "MF=%d,Nx=%g,Ny=%g,Nz=%g,Ux=%g,Uy=%g,Uz=%g",
433     curparams.hsiz,
434     curparams.nrm[0], curparams.nrm[1], curparams.nrm[2],
435     curparams.vup[0], curparams.vup[1], curparams.vup[2]);
436     params = savqstr(sbuf);
437     binv = "rbin";
438     nbins = "Nrbins";
439 greg 2.1 } else if (tolower(curparams.hemis[0]) == 'k' &&
440     !curparams.hemis[1] |
441     (tolower(curparams.hemis[1]) == 'f') |
442     (curparams.hemis[1] == '1')) {
443 greg 2.2 calfn = kfullfn; kfullfn = NULL;
444     binf = "kbin";
445     nbins = "Nkbins";
446 greg 2.1 } else if (tolower(curparams.hemis[0]) == 'k' &&
447     (tolower(curparams.hemis[1]) == 'h') |
448     (curparams.hemis[1] == '2')) {
449 greg 2.2 calfn = khalffn; khalffn = NULL;
450     binf = "khbin";
451     nbins = "Nkhbins";
452 greg 2.1 } else if (tolower(curparams.hemis[0]) == 'k' &&
453     (tolower(curparams.hemis[1]) == 'q') |
454     (curparams.hemis[1] == '4')) {
455 greg 2.2 calfn = kquarterfn; kquarterfn = NULL;
456     binf = "kqbin";
457     nbins = "Nkqbins";
458 greg 2.1 } else {
459     fprintf(stderr, "%s: unrecognized hemisphere sampling: h=%s\n",
460     progname, curparams.hemis);
461     exit(1);
462     }
463 greg 2.4 if (!uniform & (curparams.slist->styp == ST_SOURCE)) {
464     SURF *sp;
465     for (sp = curparams.slist; sp != NULL; sp = sp->next)
466     if (fabs(sp->area - PI) > 1e-3) {
467     fprintf(stderr, "%s: source '%s' must be 180-degrees\n",
468     progname, sp->sname);
469     exit(1);
470     }
471     }
472 greg 2.2 if (calfn != NULL) { /* add cal file if needed */
473     CHECKARGC(2);
474     rcarg[nrcargs++] = "-f";
475     rcarg[nrcargs++] = calfn;
476     }
477 greg 2.3 if (params != NULL) { /* parameters _after_ cal file */
478     CHECKARGC(2);
479     rcarg[nrcargs++] = "-p";
480     rcarg[nrcargs++] = params;
481     }
482 greg 2.2 if (nbins != NULL) { /* add #bins if set */
483     CHECKARGC(2);
484     rcarg[nrcargs++] = "-bn";
485     rcarg[nrcargs++] = nbins;
486     }
487 greg 2.3 if (binv != NULL) {
488 greg 2.2 CHECKARGC(2); /* assign bin variable */
489     rcarg[nrcargs++] = "-b";
490     rcarg[nrcargs++] = binv;
491     } else if (binf != NULL) {
492     CHECKARGC(2); /* assign bin function */
493 greg 2.3 rcarg[nrcargs++] = "-b";
494 greg 2.2 sprintf(sbuf, "%s(%g,%g,%g,%g,%g,%g)", binf,
495     curparams.nrm[0], curparams.nrm[1], curparams.nrm[2],
496     curparams.vup[0], curparams.vup[1], curparams.vup[2]);
497     rcarg[nrcargs++] = savqstr(sbuf);
498     }
499     CHECKARGC(2); /* modifier argument goes last */
500 greg 2.1 rcarg[nrcargs++] = "-m";
501     rcarg[nrcargs++] = savqstr(curmod);
502     }
503    
504     /* Make randomly oriented tangent plane axes for given normal direction */
505     static void
506     make_axes(FVECT uva[2], const FVECT nrm)
507     {
508     int i;
509    
510     uva[1][0] = 0.5 - frandom();
511     uva[1][1] = 0.5 - frandom();
512     uva[1][2] = 0.5 - frandom();
513     for (i = 3; i--; )
514     if ((-0.6 < nrm[i]) & (nrm[i] < 0.6))
515     break;
516     if (i < 0) {
517     fputs(progname, stderr);
518     fputs(": bad surface normal in make_axes!\n", stderr);
519     exit(1);
520     }
521     uva[1][i] = 1.0;
522     VCROSS(uva[0], uva[1], nrm);
523     normalize(uva[0]);
524     VCROSS(uva[1], nrm, uva[0]);
525     }
526    
527     /* Illegal sender surfaces end up here */
528     static int
529     ssamp_bad(FVECT orig, SURF *sp, double x)
530     {
531     fprintf(stderr, "%s: illegal sender surface '%s'\n",
532     progname, sp->sname);
533     return(0);
534     }
535    
536     /* Generate origin on ring surface from uniform random variable */
537     static int
538     ssamp_ring(FVECT orig, SURF *sp, double x)
539     {
540     FVECT *uva = (FVECT *)sp->priv;
541     double samp2[2];
542     double uv[2];
543     int i;
544    
545     if (uva == NULL) { /* need tangent axes */
546     uva = (FVECT *)malloc(sizeof(FVECT)*2);
547     if (uva == NULL) {
548     fputs(progname, stderr);
549     fputs(": out of memory in ssamp_ring!\n", stderr);
550     return(0);
551     }
552     make_axes(uva, sp->snrm);
553     sp->priv = (void *)uva;
554     }
555     SDmultiSamp(samp2, 2, x);
556 greg 2.6 samp2[0] = sqrt(samp2[0]*sp->area*(1./PI) + sp->farg[6]*sp->farg[6]);
557 greg 2.1 samp2[1] *= 2.*PI;
558     uv[0] = samp2[0]*tcos(samp2[1]);
559     uv[1] = samp2[0]*tsin(samp2[1]);
560     for (i = 3; i--; )
561     orig[i] = sp->farg[i] + uv[0]*uva[0][i] + uv[1]*uva[1][i];
562     return(1);
563     }
564    
565     /* Add triangle to polygon's list (call-back function) */
566     static int
567     add_triangle(const Vert2_list *tp, int a, int b, int c)
568     {
569     POLYTRIS *ptp = (POLYTRIS *)tp->p;
570     struct ptri *trip = ptp->tri + ptp->ntris++;
571    
572     trip->vndx[0] = a;
573     trip->vndx[1] = b;
574     trip->vndx[2] = c;
575     return(1);
576     }
577    
578     /* Generate origin on polygon surface from uniform random variable */
579     static int
580     ssamp_poly(FVECT orig, SURF *sp, double x)
581     {
582     POLYTRIS *ptp = (POLYTRIS *)sp->priv;
583     double samp2[2];
584     double *v0, *v1, *v2;
585     int i;
586    
587     if (ptp == NULL) { /* need to triangulate */
588     ptp = (POLYTRIS *)malloc(sizeof(POLYTRIS) +
589     sizeof(struct ptri)*(sp->nfargs/3 - 3));
590     if (ptp == NULL)
591     goto memerr;
592     if (sp->nfargs == 3) { /* simple case */
593     ptp->ntris = 1;
594     ptp->tri[0].vndx[0] = 0;
595     ptp->tri[0].vndx[1] = 1;
596     ptp->tri[0].vndx[2] = 2;
597     ptp->tri[0].afrac = 1;
598     } else {
599     Vert2_list *v2l = polyAlloc(sp->nfargs/3);
600     if (v2l == NULL)
601     goto memerr;
602     make_axes(ptp->uva, sp->snrm);
603     for (i = v2l->nv; i--; ) {
604     v2l->v[i].mX = DOT(sp->farg+3*i, ptp->uva[0]);
605     v2l->v[i].mY = DOT(sp->farg+3*i, ptp->uva[1]);
606     }
607     ptp->ntris = 0;
608     v2l->p = (void *)ptp;
609 greg 2.7 if (!polyTriangulate(v2l, add_triangle)) {
610     fprintf(stderr,
611     "%s: cannot triangulate polygon '%s'\n",
612     progname, sp->sname);
613 greg 2.1 return(0);
614 greg 2.7 }
615 greg 2.1 for (i = ptp->ntris; i--; ) {
616     int a = ptp->tri[i].vndx[0];
617     int b = ptp->tri[i].vndx[1];
618     int c = ptp->tri[i].vndx[2];
619     ptp->tri[i].afrac =
620     (v2l->v[a].mX*v2l->v[b].mY -
621     v2l->v[b].mX*v2l->v[a].mY +
622     v2l->v[b].mX*v2l->v[c].mY -
623     v2l->v[c].mX*v2l->v[b].mY +
624     v2l->v[c].mX*v2l->v[a].mY -
625     v2l->v[a].mX*v2l->v[c].mY) /
626     (2.*sp->area);
627     }
628     polyFree(v2l);
629     }
630     sp->priv = (void *)ptp;
631     }
632     /* pick triangle by partial area */
633     for (i = 0; i < ptp->ntris && x > ptp->tri[i].afrac; i++)
634     x -= ptp->tri[i].afrac;
635     SDmultiSamp(samp2, 2, x/ptp->tri[i].afrac);
636     samp2[0] *= samp2[1] = sqrt(samp2[1]);
637     samp2[1] = 1. - samp2[1];
638     v0 = sp->farg + 3*ptp->tri[i].vndx[0];
639     v1 = sp->farg + 3*ptp->tri[i].vndx[1];
640     v2 = sp->farg + 3*ptp->tri[i].vndx[2];
641     for (i = 3; i--; )
642     orig[i] = v0[i] + samp2[0]*(v1[i] - v0[i])
643     + samp2[1]*(v2[i] - v0[i]) ;
644     return(1);
645     memerr:
646     fputs(progname, stderr);
647     fputs(": out of memory in ssamp_poly!\n", stderr);
648     return(0);
649     }
650    
651     /* Compute sample origin based on projected areas of sender subsurfaces */
652     static int
653     sample_origin(PARAMS *p, FVECT orig, const FVECT rdir, double x)
654     {
655     static double *projsa;
656     static int nall;
657     double tarea = 0;
658     int i;
659     SURF *sp;
660     /* special case for lone surface */
661     if (p->nsurfs == 1) {
662     sp = p->slist;
663 greg 2.7 if (DOT(sp->snrm, rdir) >= -FTINY) {
664     fprintf(stderr,
665     "%s: internal - sample behind sender '%s'\n",
666     progname, sp->sname);
667     return(0);
668     }
669 greg 2.1 return((*orig_in_surf[sp->styp])(orig, sp, x));
670     }
671     if (p->nsurfs > nall) { /* (re)allocate surface area cache */
672     if (projsa) free(projsa);
673     projsa = (double *)malloc(sizeof(double)*p->nsurfs);
674     if (!projsa) return(0);
675     nall = p->nsurfs;
676     }
677     /* compute projected areas */
678     for (i = 0, sp = p->slist; sp != NULL; i++, sp = sp->next) {
679     projsa[i] = -DOT(sp->snrm, rdir) * sp->area;
680     tarea += projsa[i] *= (double)(projsa[i] > FTINY);
681     }
682 greg 2.7 if (tarea <= FTINY) { /* wrong side of sender? */
683     fputs(progname, stderr);
684     fputs(": internal - sample behind all sender elements!\n",
685     stderr);
686 greg 2.1 return(0);
687 greg 2.7 }
688 greg 2.1 tarea *= x; /* get surface from list */
689     for (i = 0, sp = p->slist; tarea > projsa[i]; sp = sp->next)
690     tarea -= projsa[i++];
691     return((*orig_in_surf[sp->styp])(orig, sp, tarea/projsa[i]));
692     }
693    
694     /* Uniform sample generator */
695     static int
696     sample_uniform(PARAMS *p, int b, FILE *fp)
697     {
698     int n = sampcnt;
699     double samp3[3];
700     double duvw[3];
701     FVECT orig_dir[2];
702     int i;
703    
704     if (fp == NULL) /* just requesting number of bins? */
705     return(1);
706    
707     while (n--) { /* stratified hemisphere sampling */
708     SDmultiSamp(samp3, 3, (n+frandom())/sampcnt);
709     SDsquare2disk(duvw, samp3[1], samp3[2]);
710     duvw[2] = -sqrt(1. - duvw[0]*duvw[0] - duvw[1]*duvw[1]);
711     for (i = 3; i--; )
712     orig_dir[1][i] = duvw[0]*p->udir[i] +
713     duvw[1]*p->vdir[i] +
714     duvw[2]*p->nrm[i] ;
715     if (!sample_origin(p, orig_dir[0], orig_dir[1], samp3[0]))
716     return(0);
717     if (fwrite(orig_dir, sizeof(FVECT), 2, fp) != 2)
718     return(0);
719     }
720     return(1);
721     }
722    
723     /* Shirly-Chiu sample generator */
724     static int
725     sample_shirchiu(PARAMS *p, int b, FILE *fp)
726     {
727     int n = sampcnt;
728     double samp3[3];
729     double duvw[3];
730     FVECT orig_dir[2];
731     int i;
732    
733     if (fp == NULL) /* just requesting number of bins? */
734     return(p->hsiz*p->hsiz);
735    
736     while (n--) { /* stratified sampling */
737     SDmultiSamp(samp3, 3, (n+frandom())/sampcnt);
738     SDsquare2disk(duvw, (b/p->hsiz + samp3[1])/curparams.hsiz,
739     (b%p->hsiz + samp3[2])/curparams.hsiz);
740     duvw[2] = sqrt(1. - duvw[0]*duvw[0] - duvw[1]*duvw[1]);
741     for (i = 3; i--; )
742     orig_dir[1][i] = -duvw[0]*p->udir[i] -
743     duvw[1]*p->vdir[i] -
744     duvw[2]*p->nrm[i] ;
745     if (!sample_origin(p, orig_dir[0], orig_dir[1], samp3[0]))
746     return(0);
747     if (fwrite(orig_dir, sizeof(FVECT), 2, fp) != 2)
748     return(0);
749     }
750     return(1);
751     }
752    
753     /* Reinhart/Tregenza sample generator */
754     static int
755     sample_reinhart(PARAMS *p, int b, FILE *fp)
756     {
757     #define T_NALT 7
758     static const int tnaz[T_NALT] = {30, 30, 24, 24, 18, 12, 6};
759     const int RowMax = T_NALT*p->hsiz + 1;
760 greg 2.7 const double RAH = (.5*PI)/(RowMax-.5);
761 greg 2.1 #define rnaz(r) (r >= RowMax-1 ? 1 : p->hsiz*tnaz[r/p->hsiz])
762     int n = sampcnt;
763     int row, col;
764     double samp3[3];
765     double alt, azi;
766     double duvw[3];
767     FVECT orig_dir[2];
768     int i;
769    
770     if (fp == NULL) { /* just requesting number of bins? */
771     n = 0;
772     for (row = RowMax; row--; ) n += rnaz(row);
773     return(n);
774     }
775     row = 0; /* identify row & column */
776     col = b;
777     while (col >= rnaz(row)) {
778     col -= rnaz(row);
779     ++row;
780     }
781     while (n--) { /* stratified sampling */
782     SDmultiSamp(samp3, 3, (n+frandom())/sampcnt);
783     alt = (row+samp3[1])*RAH;
784     azi = (2.*PI)*(col+samp3[2]-.5)/rnaz(row);
785 greg 2.7 duvw[2] = cos(alt); /* measured from horizon */
786 greg 2.1 duvw[0] = tcos(azi)*duvw[2];
787     duvw[1] = tsin(azi)*duvw[2];
788     duvw[2] = sqrt(1. - duvw[2]*duvw[2]);
789     for (i = 3; i--; )
790     orig_dir[1][i] = -duvw[0]*p->udir[i] -
791     duvw[1]*p->vdir[i] -
792     duvw[2]*p->nrm[i] ;
793     if (!sample_origin(p, orig_dir[0], orig_dir[1], samp3[0]))
794     return(0);
795     if (fwrite(orig_dir, sizeof(FVECT), 2, fp) != 2)
796     return(0);
797     }
798     return(1);
799     #undef rnaz
800     #undef T_NALT
801     }
802    
803     /* Klems sample generator */
804     static int
805     sample_klems(PARAMS *p, int b, FILE *fp)
806     {
807     static const char bname[4][20] = {
808     "LBNL/Klems Full",
809     "LBNL/Klems Half",
810     "INTERNAL ERROR",
811     "LBNL/Klems Quarter"
812     };
813     static ANGLE_BASIS *kbasis[4];
814     const int bi = p->hemis[1] - '1';
815     int n = sampcnt;
816     double samp2[2];
817     double duvw[3];
818     FVECT orig_dir[2];
819     int i;
820    
821     if (!kbasis[bi]) { /* need to get basis, first */
822     for (i = 4; i--; )
823     if (!strcasecmp(abase_list[i].name, bname[bi])) {
824     kbasis[bi] = &abase_list[i];
825     break;
826     }
827     if (i < 0) {
828     fprintf(stderr, "%s: unknown hemisphere basis '%s'\n",
829     progname, bname[bi]);
830     return(0);
831     }
832     }
833     if (fp == NULL) /* just requesting number of bins? */
834     return(kbasis[bi]->nangles);
835    
836     while (n--) { /* stratified sampling */
837     SDmultiSamp(samp2, 2, (n+frandom())/sampcnt);
838     if (!bo_getvec(duvw, b+samp2[1], kbasis[bi]))
839     return(0);
840     for (i = 3; i--; )
841     orig_dir[1][i] = duvw[0]*p->udir[i] +
842     duvw[1]*p->vdir[i] +
843     duvw[2]*p->nrm[i] ;
844     if (!sample_origin(p, orig_dir[0], orig_dir[1], samp2[0]))
845     return(0);
846     if (fwrite(orig_dir, sizeof(FVECT), 2, fp) != 2)
847     return(0);
848     }
849     return(1);
850     }
851    
852     /* Prepare hemisphere basis sampler that will send rays to rcontrib */
853     static int
854     prepare_sampler(void)
855     {
856     if (curparams.slist == NULL) { /* missing sample surface! */
857     fputs(progname, stderr);
858     fputs(": no sender surface!\n", stderr);
859     return(-1);
860     }
861     if (curparams.outfn != NULL) /* misplaced output file spec. */
862     fprintf(stderr, "%s: warning - ignoring output file in sender ('%s')\n",
863     progname, curparams.outfn);
864     /* check/set basis hemisphere */
865     if (!curparams.hemis[0]) {
866     fputs(progname, stderr);
867     fputs(": missing sender sampling type!\n", stderr);
868     return(-1);
869     }
870     if (normalize(curparams.nrm) == 0) {
871     fputs(progname, stderr);
872     fputs(": undefined normal for sender sampling\n", stderr);
873     return(-1);
874     }
875     if (normalize(curparams.vup) == 0)
876     if (fabs(curparams.nrm[2]) < .7)
877     curparams.vup[2] = 1;
878     else
879     curparams.vup[1] = 1;
880     VCROSS(curparams.udir, curparams.vup, curparams.nrm);
881     if (normalize(curparams.udir) == 0) {
882     fputs(progname, stderr);
883     fputs(": up vector coincides with sender normal\n", stderr);
884     return(-1);
885     }
886     VCROSS(curparams.vdir, curparams.nrm, curparams.udir);
887     if (tolower(curparams.hemis[0]) == 'u' | curparams.hemis[0] == '1')
888     curparams.sample_basis = sample_uniform;
889     else if (tolower(curparams.hemis[0]) == 's' &&
890     tolower(curparams.hemis[1]) == 'c')
891     curparams.sample_basis = sample_shirchiu;
892     else if ((tolower(curparams.hemis[0]) == 'r') |
893     (tolower(curparams.hemis[0]) == 't'))
894     curparams.sample_basis = sample_reinhart;
895     else if (tolower(curparams.hemis[0]) == 'k') {
896     switch (curparams.hemis[1]) {
897     case '1':
898     case '2':
899     case '4':
900     break;
901     case 'f':
902     case 'F':
903     case '\0':
904     curparams.hemis[1] = '1';
905     break;
906     case 'h':
907     case 'H':
908     curparams.hemis[1] = '2';
909     break;
910     case 'q':
911     case 'Q':
912     curparams.hemis[1] = '4';
913     break;
914     default:
915     goto unrecognized;
916     }
917     curparams.hemis[2] = '\0';
918     curparams.sample_basis = sample_klems;
919     } else
920     goto unrecognized;
921     /* return number of bins */
922     return((*curparams.sample_basis)(&curparams,0,NULL));
923     unrecognized:
924     fprintf(stderr, "%s: unrecognized sender sampling: h=%s\n",
925     progname, curparams.hemis);
926     return(-1);
927     }
928    
929     /* Compute normal and area for polygon */
930     static int
931     finish_polygon(SURF *p)
932     {
933     const int nv = p->nfargs / 3;
934     FVECT e1, e2, vc;
935     int i;
936    
937     memset(p->snrm, 0, sizeof(FVECT));
938     VSUB(e1, p->farg+3, p->farg);
939     for (i = 2; i < nv; i++) {
940     VSUB(e2, p->farg+3*i, p->farg);
941     VCROSS(vc, e1, e2);
942     p->snrm[0] += vc[0];
943     p->snrm[1] += vc[1];
944     p->snrm[2] += vc[2];
945     VCOPY(e1, e2);
946     }
947     p->area = normalize(p->snrm)*0.5;
948     return(p->area > FTINY);
949     }
950    
951     /* Add a surface to our current parameters */
952     static void
953     add_surface(int st, const char *oname, FILE *fp)
954     {
955     SURF *snew;
956     int n;
957     /* get floating-point arguments */
958     if (!fscanf(fp, "%d", &n)) return;
959     while (n-- > 0) fscanf(fp, "%*s");
960     if (!fscanf(fp, "%d", &n)) return;
961     while (n-- > 0) fscanf(fp, "%*d");
962     if (!fscanf(fp, "%d", &n) || n <= 0) return;
963     snew = (SURF *)malloc(sizeof(SURF) + sizeof(double)*(n-1));
964     if (snew == NULL) {
965     fputs(progname, stderr);
966     fputs(": out of memory!\n", stderr);
967     exit(1);
968     }
969     strncpy(snew->sname, oname, sizeof(snew->sname)-1);
970     snew->sname[sizeof(snew->sname)-1] = '\0';
971     snew->styp = st;
972     snew->priv = NULL;
973     snew->nfargs = n;
974     for (n = 0; n < snew->nfargs; n++)
975     if (fscanf(fp, "%lf", &snew->farg[n]) != 1) {
976     fprintf(stderr, "%s: error reading arguments for '%s'\n",
977     progname, oname);
978     exit(1);
979     }
980     switch (st) {
981     case ST_RING:
982     if (snew->nfargs != 8)
983     goto badcount;
984     VCOPY(snew->snrm, snew->farg+3);
985     if (normalize(snew->snrm) == 0)
986     goto badnorm;
987     if (snew->farg[7] < snew->farg[6]) {
988     double t = snew->farg[7];
989     snew->farg[7] = snew->farg[6];
990     snew->farg[6] = t;
991     }
992     snew->area = PI*(snew->farg[7]*snew->farg[7] -
993     snew->farg[6]*snew->farg[6]);
994     break;
995     case ST_POLY:
996     if (snew->nfargs < 9 || snew->nfargs % 3)
997     goto badcount;
998     finish_polygon(snew);
999     break;
1000     case ST_SOURCE:
1001     if (snew->nfargs != 4)
1002     goto badcount;
1003 greg 2.8 for (n = 3; n--; ) /* need to reverse "normal" */
1004     snew->snrm[n] = -snew->farg[n];
1005 greg 2.1 if (normalize(snew->snrm) == 0)
1006     goto badnorm;
1007 greg 2.4 snew->area = sin((PI/180./2.)*snew->farg[3]);
1008     snew->area *= PI*snew->area;
1009 greg 2.1 break;
1010     }
1011     if (snew->area <= FTINY) {
1012     fprintf(stderr, "%s: warning - zero area for surface '%s'\n",
1013     progname, oname);
1014     free(snew);
1015     return;
1016     }
1017     VSUM(curparams.nrm, curparams.nrm, snew->snrm, snew->area);
1018     snew->next = curparams.slist;
1019     curparams.slist = snew;
1020     curparams.nsurfs++;
1021     return;
1022     badcount:
1023 greg 2.7 fprintf(stderr, "%s: bad argument count for surface element '%s'\n",
1024 greg 2.1 progname, oname);
1025     exit(1);
1026     badnorm:
1027 greg 2.7 fprintf(stderr, "%s: bad orientation for surface element '%s'\n",
1028 greg 2.1 progname, oname);
1029     exit(1);
1030     }
1031    
1032     /* Parse a receiver object (look for modifiers to add to rcontrib command) */
1033     static int
1034     add_recv_object(FILE *fp)
1035     {
1036     int st;
1037     char thismod[128], otype[32], oname[128];
1038     int n;
1039    
1040     if (fscanf(fp, "%s %s %s", thismod, otype, oname) != 3)
1041     return(0); /* must have hit EOF! */
1042     if (!strcmp(otype, "alias")) {
1043     fscanf(fp, "%*s"); /* skip alias */
1044     return(0);
1045     }
1046     /* is it a new receiver? */
1047     if ((st = surf_type(otype)) != ST_NONE) {
1048     if (curparams.slist != NULL && (st == ST_SOURCE) ^
1049     (curparams.slist->styp == ST_SOURCE)) {
1050     fputs(progname, stderr);
1051     fputs(": cannot mix source/non-source receivers!\n", stderr);
1052     return(-1);
1053     }
1054     if (strcmp(thismod, curmod)) {
1055     if (curmod[0]) { /* output last receiver? */
1056     finish_receiver();
1057     clear_params(&curparams, 1);
1058     }
1059 greg 2.6 parse_params(&curparams, newparams);
1060     newparams[0] = '\0';
1061 greg 2.1 strcpy(curmod, thismod);
1062     }
1063     add_surface(st, oname, fp); /* read & store surface */
1064     return(1);
1065     }
1066     /* else skip arguments */
1067 greg 2.8 if (!fscanf(fp, "%d", &n)) return(0);
1068 greg 2.1 while (n-- > 0) fscanf(fp, "%*s");
1069     if (!fscanf(fp, "%d", &n)) return;
1070     while (n-- > 0) fscanf(fp, "%*d");
1071     if (!fscanf(fp, "%d", &n)) return;
1072     while (n-- > 0) fscanf(fp, "%*f");
1073     return(0);
1074     }
1075    
1076     /* Parse a sender object */
1077     static int
1078     add_send_object(FILE *fp)
1079     {
1080     int st;
1081     char otype[32], oname[128];
1082     int n;
1083    
1084     if (fscanf(fp, "%*s %s %s", otype, oname) != 2)
1085     return(0); /* must have hit EOF! */
1086     if (!strcmp(otype, "alias")) {
1087     fscanf(fp, "%*s"); /* skip alias */
1088     return(0);
1089     }
1090     /* is it a new surface? */
1091     if ((st = surf_type(otype)) != ST_NONE) {
1092     if (st == ST_SOURCE) {
1093     fputs(progname, stderr);
1094     fputs(": cannot use source as a sender!\n", stderr);
1095     return(-1);
1096     }
1097 greg 2.6 parse_params(&curparams, newparams);
1098     newparams[0] = '\0';
1099 greg 2.1 add_surface(st, oname, fp); /* read & store surface */
1100     return(0);
1101     }
1102     /* else skip arguments */
1103 greg 2.8 if (!fscanf(fp, "%d", &n)) return(0);
1104 greg 2.1 while (n-- > 0) fscanf(fp, "%*s");
1105     if (!fscanf(fp, "%d", &n)) return;
1106     while (n-- > 0) fscanf(fp, "%*d");
1107     if (!fscanf(fp, "%d", &n)) return;
1108     while (n-- > 0) fscanf(fp, "%*f");
1109     return(0);
1110     }
1111    
1112     /* Load a Radiance scene using the given callback function for objects */
1113     static int
1114     load_scene(const char *inspec, int (*ocb)(FILE *))
1115     {
1116     int rv = 0;
1117     char inpbuf[1024];
1118     FILE *fp;
1119     int c;
1120    
1121     if (*inspec == '!')
1122     fp = popen(inspec+1, "r");
1123     else
1124     fp = fopen(inspec, "r");
1125     if (fp == NULL) {
1126     fprintf(stderr, "%s: cannot load '%s'\n", progname, inspec);
1127     return(-1);
1128     }
1129     while ((c = getc(fp)) != EOF) { /* load receiver data */
1130     if (isspace(c)) /* skip leading white space */
1131     continue;
1132     if (c == '!') { /* read from a new command */
1133     inpbuf[0] = c;
1134     if (fgetline(inpbuf+1, sizeof(inpbuf)-1, fp) != NULL) {
1135     if ((c = load_scene(inpbuf, ocb)) < 0)
1136     return(c);
1137     rv += c;
1138     }
1139     continue;
1140     }
1141     if (c == '#') { /* parameters/comment */
1142 greg 2.5 if ((c = getc(fp)) == EOF || ungetc(c,fp) == EOF)
1143     break;
1144     if (!isspace(c) && fscanf(fp, "%s", inpbuf) == 1 &&
1145 greg 2.1 !strcmp(inpbuf, PARAMSTART)) {
1146     if (fgets(inpbuf, sizeof(inpbuf), fp) != NULL)
1147 greg 2.6 strcat(newparams, inpbuf);
1148 greg 2.1 continue;
1149     }
1150     while ((c = getc(fp)) != EOF && c != '\n');
1151     ; /* else skipping comment */
1152     continue;
1153     }
1154     ungetc(c, fp); /* else check object for receiver */
1155     c = (*ocb)(fp);
1156     if (c < 0)
1157     return(c);
1158     rv += c;
1159     }
1160     /* close our input stream */
1161     c = (*inspec == '!') ? pclose(fp) : fclose(fp);
1162     if (c != 0) {
1163     fprintf(stderr, "%s: error loading '%s'\n", progname, inspec);
1164     return(-1);
1165     }
1166     return(rv);
1167     }
1168    
1169     /* Get command arguments and run program */
1170     int
1171     main(int argc, char *argv[])
1172     {
1173     char fmtopt[6] = "-faa"; /* default output is ASCII */
1174 greg 2.4 char *xrs=NULL, *yrs=NULL, *ldopt=NULL;
1175 greg 2.8 int wantIrradiance = 0;
1176 greg 2.1 char *sendfn;
1177     char sampcntbuf[32], nsbinbuf[32];
1178     FILE *rcfp;
1179     int nsbins;
1180     int a, i;
1181     /* screen rcontrib options */
1182     progname = argv[0];
1183     for (a = 1; a < argc-2 && argv[a][0] == '-'; a++) {
1184     int na = 1; /* !! Keep consistent !! */
1185     switch (argv[a][1]) {
1186     case 'v': /* verbose mode */
1187     verbose = !verbose;
1188     na = 0;
1189     continue;
1190     case 'f': /* special case for -fo, -ff, etc. */
1191     switch (argv[a][2]) {
1192     case '\0': /* cal file */
1193     goto userr;
1194     case 'o': /* force output */
1195     goto userr;
1196     case 'a': /* output format */
1197     case 'f':
1198     case 'd':
1199     case 'c':
1200 greg 2.6 if (!(fmtopt[3] = argv[a][3]))
1201     fmtopt[3] = argv[a][2];
1202     fmtopt[2] = argv[a][2];
1203 greg 2.1 na = 0;
1204     continue; /* will pass later */
1205     default:
1206     goto userr;
1207     }
1208     break;
1209 greg 2.4 case 'x': /* x-resolution */
1210     xrs = argv[++a];
1211     na = 0;
1212     continue;
1213     case 'y': /* y-resolution */
1214     yrs = argv[++a];
1215     na = 0;
1216     continue;
1217 greg 2.1 case 'c': /* number of samples */
1218 greg 2.6 sampcnt = atoi(argv[++a]);
1219 greg 2.1 if (sampcnt <= 0)
1220     goto userr;
1221     na = 0; /* we re-add this later */
1222     continue;
1223 greg 2.8 case 'I': /* only for pass-through mode */
1224     wantIrradiance = 1;
1225     na = 0;
1226     continue;
1227 greg 2.1 case 'V': /* options without arguments */
1228     case 'w':
1229     case 'u':
1230     case 'i':
1231     case 'h':
1232 greg 2.4 case 'r':
1233 greg 2.1 break;
1234     case 'n': /* options with 1 argument */
1235     case 's':
1236     case 'o':
1237     na = 2;
1238     break;
1239     case 'b': /* special case */
1240     if (argv[a][2] != 'v') goto userr;
1241     break;
1242     case 'l': /* special case */
1243 greg 2.4 if (argv[a][2] == 'd') {
1244     ldopt = argv[a];
1245     na = 0;
1246     continue;
1247     }
1248 greg 2.1 na = 2;
1249     break;
1250     case 'd': /* special case */
1251     if (argv[a][2] != 'v') na = 2;
1252     break;
1253     case 'a': /* special case */
1254     na = (argv[a][2] == 'v') ? 4 : 2;
1255     break;
1256     case 'm': /* special case */
1257     if (!argv[a][2]) goto userr;
1258     na = (argv[a][2] == 'e') | (argv[a][2] == 'a') ? 4 : 2;
1259     break;
1260     case '\0': /* pass-through mode */
1261     goto done_opts;
1262     default: /* anything else is verbotten */
1263     goto userr;
1264     }
1265     if (na <= 0) continue;
1266     CHECKARGC(na); /* pass on option */
1267     rcarg[nrcargs++] = argv[a];
1268     while (--na) /* + arguments if any */
1269     rcarg[nrcargs++] = argv[++a];
1270     }
1271     done_opts:
1272     if (a > argc-2)
1273     goto userr; /* check at end of options */
1274     sendfn = argv[a++]; /* assign sender & receiver inputs */
1275     if (sendfn[0] == '-') { /* user wants pass-through mode? */
1276     if (sendfn[1]) goto userr;
1277     sendfn = NULL;
1278 greg 2.8 if (wantIrradiance) {
1279     CHECKARGC(1);
1280     rcarg[nrcargs++] = "-I";
1281     }
1282 greg 2.4 if (xrs) {
1283     CHECKARGC(2);
1284     rcarg[nrcargs++] = "-x";
1285     rcarg[nrcargs++] = xrs;
1286     }
1287     if (yrs) {
1288     CHECKARGC(2);
1289     rcarg[nrcargs++] = "-y";
1290     rcarg[nrcargs++] = yrs;
1291     }
1292     if (ldopt) {
1293     CHECKARGC(1);
1294     rcarg[nrcargs++] = ldopt;
1295     }
1296 greg 2.1 if (sampcnt <= 0) sampcnt = 1;
1297 greg 2.8 } else { /* else in sampling mode */
1298     if (wantIrradiance) {
1299     fputs(progname, stderr);
1300     fputs(": -I supported for pass-through only\n", stderr);
1301     return(1);
1302     }
1303 greg 2.6 fmtopt[2] = (sizeof(RREAL)==sizeof(double)) ? 'd' : 'f';
1304 greg 2.1 if (sampcnt <= 0) sampcnt = 10000;
1305     }
1306     sprintf(sampcntbuf, "%d", sampcnt);
1307     CHECKARGC(3); /* add our format & sample count */
1308     rcarg[nrcargs++] = fmtopt;
1309     rcarg[nrcargs++] = "-c";
1310     rcarg[nrcargs++] = sampcntbuf;
1311     /* add receiver arguments to rcontrib */
1312     if (load_scene(argv[a], add_recv_object) < 0)
1313     return(1);
1314     finish_receiver();
1315     if (sendfn == NULL) { /* pass-through mode? */
1316     CHECKARGC(1); /* add octree */
1317     rcarg[nrcargs++] = oconv_command(argc-a, argv+a);
1318     rcarg[nrcargs] = NULL;
1319     return(my_exec(rcarg)); /* rcontrib does everything */
1320     }
1321     clear_params(&curparams, 0); /* else load sender surface & params */
1322     if (load_scene(sendfn, add_send_object) < 0)
1323     return(1);
1324     if ((nsbins = prepare_sampler()) <= 0)
1325     return(1);
1326     CHECKARGC(3); /* add row count and octree */
1327     rcarg[nrcargs++] = "-y";
1328     sprintf(nsbinbuf, "%d", nsbins);
1329     rcarg[nrcargs++] = nsbinbuf;
1330     rcarg[nrcargs++] = oconv_command(argc-a, argv+a);
1331     rcarg[nrcargs] = NULL;
1332     /* open pipe to rcontrib process */
1333     if ((rcfp = popen_arglist(rcarg, "w")) == NULL)
1334     return(1);
1335     SET_FILE_BINARY(rcfp);
1336     #ifdef getc_unlocked
1337     flockfile(rcfp);
1338     #endif
1339     if (verbose) {
1340     fprintf(stderr, "%s: sampling %d directions", progname, nsbins);
1341     if (curparams.nsurfs > 1)
1342 greg 2.7 fprintf(stderr, " (%d elements)\n", curparams.nsurfs);
1343 greg 2.1 else
1344     fputc('\n', stderr);
1345     }
1346     for (i = 0; i < nsbins; i++) /* send rcontrib ray samples */
1347     if (!(*curparams.sample_basis)(&curparams, i, rcfp))
1348     return(1);
1349     return(pclose(rcfp) == 0); /* all finished! */
1350     userr:
1351     if (a < argc-2)
1352     fprintf(stderr, "%s: unsupported option '%s'", progname, argv[a]);
1353 greg 2.5 fprintf(stderr, "Usage: %s [-v][rcontrib options] sender.rad receiver.rad [system.rad ..]\n",
1354 greg 2.1 progname);
1355     return(1);
1356     }