--- ray/src/rt/ambient.c 2016/04/24 16:21:32 2.102 +++ ray/src/rt/ambient.c 2023/01/27 19:57:08 2.116 @@ -1,4 +1,4 @@ -static const char RCSid[] = "$Id: ambient.c,v 2.102 2016/04/24 16:21:32 greg Exp $"; +static const char RCSid[] = "$Id: ambient.c,v 2.116 2023/01/27 19:57:08 greg Exp $"; /* * ambient.c - routines dealing with ambient (inter-reflected) component. * @@ -12,6 +12,7 @@ static const char RCSid[] = "$Id: ambient.c,v 2.102 20 #include "platform.h" #include "ray.h" #include "otypes.h" +#include "otspecial.h" #include "resolu.h" #include "ambient.h" #include "random.h" @@ -21,8 +22,6 @@ static const char RCSid[] = "$Id: ambient.c,v 2.102 20 #define OCTSCALE 1.0 /* ceil((valid rad.)/(cube size)) */ #endif -extern char *shm_boundary; /* memory sharing boundary */ - #ifndef MAXASET #define MAXASET 4095 /* maximum number of elements in ambient set */ #endif @@ -36,47 +35,19 @@ static AMBTREE atrunk; /* our ambient trunk node */ static FILE *ambfp = NULL; /* ambient file pointer */ static int nunflshed = 0; /* number of unflushed ambient values */ -#ifndef SORT_THRESH -#ifdef SMLMEM -#define SORT_THRESH ((16L<<20)/sizeof(AMBVAL)) -#else -#define SORT_THRESH ((64L<<20)/sizeof(AMBVAL)) -#endif -#endif -#ifndef SORT_INTVL -#define SORT_INTVL (SORT_THRESH<<1) -#endif -#ifndef MAX_SORT_INTVL -#define MAX_SORT_INTVL (SORT_INTVL<<6) -#endif - - static double avsum = 0.; /* computed ambient value sum (log) */ static unsigned int navsum = 0; /* number of values in avsum */ static unsigned int nambvals = 0; /* total number of indirect values */ static unsigned int nambshare = 0; /* number of values from file */ -static unsigned long ambclock = 0; /* ambient access clock */ -static unsigned long lastsort = 0; /* time of last value sort */ -static long sortintvl = SORT_INTVL; /* time until next sort */ static FILE *ambinp = NULL; /* auxiliary file for input */ static long lastpos = -1; /* last flush position */ -#define MAXACLOCK (1L<<30) /* clock turnover value */ - /* - * Track access times unless we are sharing ambient values - * through memory on a multiprocessor, when we want to avoid - * claiming our own memory (copy on write). Go ahead anyway - * if more than two thirds of our values are unshared. - * Compile with -Dtracktime=0 to turn this code off. - */ -#ifndef tracktime -#define tracktime (shm_boundary == NULL || nambvals > 3*nambshare) -#endif - #define AMBFLUSH (BUFSIZ/AMBVALSIZ) #define newambval() (AMBVAL *)malloc(sizeof(AMBVAL)) +#define tfunc(lwr, x, upr) (((x)-(lwr))/((upr)-(lwr))) + static void initambfile(int creat); static void avsave(AMBVAL *av); static AMBVAL *avstore(AMBVAL *aval); @@ -85,14 +56,18 @@ static void freeambtree(AMBTREE *atp); typedef void unloadtf_t(AMBVAL *); static unloadtf_t avinsert; -static unloadtf_t av2list; static unloadtf_t avfree; static void unloadatree(AMBTREE *at, unloadtf_t *f); -static int aposcmp(const void *avp1, const void *avp2); -static int avlmemi(AMBVAL *avaddr); -static void sortambvals(int always); +static void sortambvals(void); +static int plugaleak(RAY *r, AMBVAL *ap, FVECT anorm, double ang); +static double sumambient(COLOR acol, RAY *r, FVECT rn, int al, + AMBTREE *at, FVECT c0, double s); +static int makeambient(COLOR acol, RAY *r, FVECT rn, int al); +static int extambient(COLOR cr, AMBVAL *ap, FVECT pv, FVECT nv, + FVECT uvw[3]); + #ifdef F_SETLKW static void aflock(int typ); #endif @@ -131,8 +106,8 @@ setambacc( /* set ambient accuracy */ newa *= (newa > 0); if (fabs(newa - olda) >= .05*(newa + olda)) { ambacc = newa; - if (nambvals > 0) - sortambvals(1); /* rebuild tree */ + if (ambacc > FTINY && nambvals > 0) + sortambvals(); /* rebuild tree */ } } @@ -219,9 +194,6 @@ ambdone(void) /* close ambient file and free memory navsum = 0; nambvals = 0; nambshare = 0; - ambclock = 0; - lastsort = 0; - sortintvl = SORT_INTVL; } @@ -254,19 +226,7 @@ ambnotify( /* record new modifier */ } } -/************ THE FOLLOWING ROUTINES DIFFER BETWEEN NEW & OLD ***************/ -#ifndef OLDAMB - -#define tfunc(lwr, x, upr) (((x)-(lwr))/((upr)-(lwr))) - -static int plugaleak(RAY *r, AMBVAL *ap, FVECT anorm, double ang); -static double sumambient(COLOR acol, RAY *r, FVECT rn, int al, - AMBTREE *at, FVECT c0, double s); -static int makeambient(COLOR acol, RAY *r, FVECT rn, int al); -static int extambient(COLOR cr, AMBVAL *ap, FVECT pv, FVECT nv, - FVECT uvw[3]); - void multambient( /* compute ambient component & multiply by coef. */ COLOR aval, @@ -274,9 +234,10 @@ multambient( /* compute ambient component & multiply FVECT nrm ) { + static double logAvgAbsorp = 1; static int rdepth = 0; /* ambient recursion */ COLOR acol, caustic; - int ok; + int i, ok; double d, l; /* PMAP: Factor in ambient from photon map, if enabled and ray is @@ -285,6 +246,9 @@ multambient( /* compute ambient component & multiply if (ambPmap(aval, r, rdepth)) return; + if (logAvgAbsorp > 0) /* exclude in -aw to avoid growth */ + logAvgAbsorp = log(1.-AVGREFL); + /* PMAP: Factor in specular-diffuse ambient (caustics) from photon * map, if enabled and ray is primary, else caustic is zero. Continue * with RADIANCE ambient calculation */ @@ -302,22 +266,33 @@ multambient( /* compute ambient component & multiply goto dumbamb; if (ambacc <= FTINY) { /* no ambient storage */ + FVECT uvd[2]; + float dgrad[2], *dgp = NULL; + + if (nrm != r->ron && DOT(nrm,r->ron) < 0.9999) + dgp = dgrad; /* compute rotational grad. */ copycolor(acol, aval); rdepth++; ok = doambient(acol, r, r->rweight, - NULL, NULL, NULL, NULL, NULL); + uvd, NULL, NULL, dgp, NULL); rdepth--; if (!ok) goto dumbamb; + if ((ok > 0) & (dgp != NULL)) { /* apply texture */ + FVECT v1; + VCROSS(v1, r->ron, nrm); + d = 1.0; + for (i = 3; i--; ) + d += v1[i] * (dgp[0]*uvd[0][i] + dgp[1]*uvd[1][i]); + if (d >= 0.05) + scalecolor(acol, d); + } copycolor(aval, acol); /* PMAP: add in caustic */ addcolor(aval, caustic); return; } - - if (tracktime) /* sort to minimize thrashing */ - sortambvals(0); /* interpolate ambient value */ setcolor(acol, 0.0, 0.0, 0.0); d = sumambient(acol, r, nrm, rdepth, @@ -356,13 +331,13 @@ dumbamb: /* return global value */ l = bright(ambval); /* average in computations */ if (l > FTINY) { - d = (log(l)*(double)ambvwt + avsum) / + d = (log(l)*(double)ambvwt + avsum + logAvgAbsorp*navsum) / (double)(ambvwt + navsum); d = exp(d) / l; scalecolor(aval, d); multcolor(aval, ambval); /* apply color of ambval */ } else { - d = exp( avsum / (double)navsum ); + d = exp( avsum/(double)navsum + logAvgAbsorp ); scalecolor(aval, d); /* neutral color */ } } @@ -401,7 +376,8 @@ plugaleak(RAY *r, AMBVAL *ap, FVECT anorm, double ang) VSUM(rtst.rdir, vdif, anorm, t[1]); /* further dist. > plane */ rtst.rmax = normalize(rtst.rdir); /* short ray test */ while (localhit(&rtst, &thescene)) { /* check for occluder */ - if (rtst.ro->omod != OVOID && + OBJREC *m = findmaterial(rtst.ro); + if (m != NULL && !istransp(m->otype) && !isBSDFproxy(m) && (rtst.clipset == NULL || !inset(rtst.clipset, rtst.ro->omod))) return(1); /* plug light leak */ @@ -459,9 +435,6 @@ sumambient( /* get interpolated ambient value */ double u, v, d, delta_r2, delta_t2; COLOR ct; FVECT uvw[3]; - /* record access */ - if (tracktime) - av->latick = ambclock; /* * Ambient level test */ @@ -483,7 +456,7 @@ sumambient( /* get interpolated ambient value */ */ VSUB(ck0, r->rop, av->pos); d = DOT(ck0, uvw[2]); - if (d < -minarad*ambacc-.001) + if (d < -minarad*ambacc) continue; d /= av->rad[0]; delta_t2 = d*d; @@ -570,6 +543,7 @@ extambient( /* extrapolate value at pv, nv */ ) { const double min_d = 0.05; + const double max_d = 20.; static FVECT my_uvw[3]; FVECT v1; int i; @@ -589,8 +563,10 @@ extambient( /* extrapolate value at pv, nv */ for (i = 3; i--; ) d += v1[i] * (ap->gdir[0]*uvw[0][i] + ap->gdir[1]*uvw[1][i]); - if (d < min_d) /* should not use if we can avoid it */ + if (d < min_d) /* clamp min/max scaling */ d = min_d; + else if (d > max_d) + d = max_d; copycolor(cr, ap->val); scalecolor(cr, d); return(d > min_d); @@ -640,337 +616,7 @@ avinsert( /* insert ambient value in our tree */ } -#else /* ! NEWAMB */ - -static double sumambient(COLOR acol, RAY *r, FVECT rn, int al, - AMBTREE *at, FVECT c0, double s); -static double makeambient(COLOR acol, RAY *r, FVECT rn, int al); -static void extambient(COLOR cr, AMBVAL *ap, FVECT pv, FVECT nv); - - -void -multambient( /* compute ambient component & multiply by coef. */ - COLOR aval, - RAY *r, - FVECT nrm -) -{ - static int rdepth = 0; /* ambient recursion */ - COLOR acol, caustic; - double d, l; - - /* PMAP: Factor in ambient from global photon map (if enabled) and return - * as all ambient components accounted for */ - if (ambPmap(aval, r, rdepth)) - return; - - /* PMAP: Otherwise factor in ambient from caustic photon map - * (ambPmapCaustic() returns zero if caustic photons disabled) and - * continue with RADIANCE ambient calculation */ - copycolor(caustic, aval); - ambPmapCaustic(caustic, r, rdepth); - - if (ambdiv <= 0) /* no ambient calculation */ - goto dumbamb; - /* check number of bounces */ - if (rdepth >= ambounce) - goto dumbamb; - /* check ambient list */ - if (ambincl != -1 && r->ro != NULL && - ambincl != inset(ambset, r->ro->omod)) - goto dumbamb; - - if (ambacc <= FTINY) { /* no ambient storage */ - copycolor(acol, aval); - rdepth++; - d = doambient(acol, r, r->rweight, NULL, NULL); - rdepth--; - if (d <= FTINY) - goto dumbamb; - copycolor(aval, acol); - - /* PMAP: add in caustic */ - addcolor(aval, caustic); - return; - } - - if (tracktime) /* sort to minimize thrashing */ - sortambvals(0); - /* interpolate ambient value */ - setcolor(acol, 0.0, 0.0, 0.0); - d = sumambient(acol, r, nrm, rdepth, - &atrunk, thescene.cuorg, thescene.cusize); - - if (d > FTINY) { - d = 1.0/d; - scalecolor(acol, d); - multcolor(aval, acol); - - /* PMAP: add in caustic */ - addcolor(aval, caustic); - return; - } - - rdepth++; /* need to cache new value */ - d = makeambient(acol, r, nrm, rdepth-1); - rdepth--; - - if (d > FTINY) { - multcolor(aval, acol); /* got new value */ - - /* PMAP: add in caustic */ - addcolor(aval, caustic); - return; - } - -dumbamb: /* return global value */ - if ((ambvwt <= 0) | (navsum == 0)) { - multcolor(aval, ambval); - - /* PMAP: add in caustic */ - addcolor(aval, caustic); - return; - } - - l = bright(ambval); /* average in computations */ - if (l > FTINY) { - d = (log(l)*(double)ambvwt + avsum) / - (double)(ambvwt + navsum); - d = exp(d) / l; - scalecolor(aval, d); - multcolor(aval, ambval); /* apply color of ambval */ - } else { - d = exp( avsum / (double)navsum ); - scalecolor(aval, d); /* neutral color */ - } -} - - -static double -sumambient( /* get interpolated ambient value */ - COLOR acol, - RAY *r, - FVECT rn, - int al, - AMBTREE *at, - FVECT c0, - double s -) -{ - double d, e1, e2, wt, wsum; - COLOR ct; - FVECT ck0; - int i; - int j; - AMBVAL *av; - - wsum = 0.0; - /* do this node */ - for (av = at->alist; av != NULL; av = av->next) { - double rn_dot = -2.0; - if (tracktime) - av->latick = ambclock; - /* - * Ambient level test. - */ - if (av->lvl > al || /* list sorted, so this works */ - (av->lvl == al) & (av->weight < 0.9*r->rweight)) - break; - /* - * Ambient radius test. - */ - VSUB(ck0, av->pos, r->rop); - e1 = DOT(ck0, ck0) / (av->rad * av->rad); - if (e1 > ambacc*ambacc*1.21) - continue; - /* - * Direction test using closest normal. - */ - d = DOT(av->dir, r->ron); - if (rn != r->ron) { - rn_dot = DOT(av->dir, rn); - if (rn_dot > 1.0-FTINY) - rn_dot = 1.0-FTINY; - if (rn_dot >= d-FTINY) { - d = rn_dot; - rn_dot = -2.0; - } - } - e2 = (1.0 - d) * r->rweight; - if (e2 < 0.0) - e2 = 0.0; - else if (e1 + e2 > ambacc*ambacc*1.21) - continue; - /* - * Ray behind test. - */ - d = 0.0; - for (j = 0; j < 3; j++) - d += (r->rop[j] - av->pos[j]) * - (av->dir[j] + r->ron[j]); - if (d*0.5 < -minarad*ambacc-.001) - continue; - /* - * Jittering final test reduces image artifacts. - */ - e1 = sqrt(e1); - e2 = sqrt(e2); - wt = e1 + e2; - if (wt > ambacc*(.9+.2*urand(9015+samplendx))) - continue; - /* - * Recompute directional error using perturbed normal - */ - if (rn_dot > 0.0) { - e2 = sqrt((1.0 - rn_dot)*r->rweight); - wt = e1 + e2; - } - if (wt <= 1e-3) - wt = 1e3; - else - wt = 1.0 / wt; - wsum += wt; - extambient(ct, av, r->rop, rn); - scalecolor(ct, wt); - addcolor(acol, ct); - } - if (at->kid == NULL) - return(wsum); - /* do children */ - s *= 0.5; - for (i = 0; i < 8; i++) { - for (j = 0; j < 3; j++) { - ck0[j] = c0[j]; - if (1<rop[j] < ck0[j] - OCTSCALE*s) - break; - if (r->rop[j] > ck0[j] + (1.0+OCTSCALE)*s) - break; - } - if (j == 3) - wsum += sumambient(acol, r, rn, al, - at->kid+i, ck0, s); - } - return(wsum); -} - - -static double -makeambient( /* make a new ambient value for storage */ - COLOR acol, - RAY *r, - FVECT rn, - int al -) -{ - AMBVAL amb; - FVECT gp, gd; - int i; - - amb.weight = 1.0; /* compute weight */ - for (i = al; i-- > 0; ) - amb.weight *= AVGREFL; - if (r->rweight < 0.1*amb.weight) /* heuristic override */ - amb.weight = 1.25*r->rweight; - setcolor(acol, AVGREFL, AVGREFL, AVGREFL); - /* compute ambient */ - amb.rad = doambient(acol, r, amb.weight, gp, gd); - if (amb.rad <= FTINY) { - setcolor(acol, 0.0, 0.0, 0.0); - return(0.0); - } - scalecolor(acol, 1./AVGREFL); /* undo assumed reflectance */ - /* store value */ - VCOPY(amb.pos, r->rop); - VCOPY(amb.dir, r->ron); - amb.lvl = al; - copycolor(amb.val, acol); - VCOPY(amb.gpos, gp); - VCOPY(amb.gdir, gd); - /* insert into tree */ - avsave(&amb); /* and save to file */ - if (rn != r->ron) - extambient(acol, &amb, r->rop, rn); /* texture */ - return(amb.rad); -} - - static void -extambient( /* extrapolate value at pv, nv */ - COLOR cr, - AMBVAL *ap, - FVECT pv, - FVECT nv -) -{ - FVECT v1; - int i; - double d; - - d = 1.0; /* zeroeth order */ - /* gradient due to translation */ - for (i = 0; i < 3; i++) - d += ap->gpos[i]*(pv[i]-ap->pos[i]); - /* gradient due to rotation */ - VCROSS(v1, ap->dir, nv); - d += DOT(ap->gdir, v1); - if (d <= 0.0) { - setcolor(cr, 0.0, 0.0, 0.0); - return; - } - copycolor(cr, ap->val); - scalecolor(cr, d); -} - - -static void -avinsert( /* insert ambient value in our tree */ - AMBVAL *av -) -{ - AMBTREE *at; - AMBVAL *ap; - AMBVAL avh; - FVECT ck0; - double s; - int branch; - int i; - - if (av->rad <= FTINY) - error(CONSISTENCY, "zero ambient radius in avinsert"); - at = &atrunk; - VCOPY(ck0, thescene.cuorg); - s = thescene.cusize; - while (s*(OCTSCALE/2) > av->rad*ambacc) { - if (at->kid == NULL) - if ((at->kid = newambtree()) == NULL) - error(SYSTEM, "out of memory in avinsert"); - s *= 0.5; - branch = 0; - for (i = 0; i < 3; i++) - if (av->pos[i] > ck0[i] + s) { - ck0[i] += s; - branch |= 1 << i; - } - at = at->kid + branch; - } - avh.next = at->alist; /* order by increasing level */ - for (ap = &avh; ap->next != NULL; ap = ap->next) - if ( ap->next->lvl > av->lvl || - (ap->next->lvl == av->lvl) & - (ap->next->weight <= av->weight) ) - break; - av->next = ap->next; - ap->next = (AMBVAL*)av; - at->alist = avh.next; -} - -#endif /* ! NEWAMB */ - -/************* FOLLOWING ROUTINES SAME FOR NEW & OLD METHODS ***************/ - -static void initambfile( /* initialize ambient file */ int cre8 ) @@ -993,6 +639,10 @@ initambfile( /* initialize ambient file */ ambvwt, ambounce, ambacc); fprintf(ambfp, "-ad %d -as %d -ar %d ", ambdiv, ambssamp, ambres); + fprintf(ambfp, "-dr %d -ds %g -dt %g -dc %g ", directrelay, + srcsizerat, shadthresh, shadcert); + fprintf(ambfp, "-ss %g -st %g -lr %d -lw %g ", specjitter, + specthresh, maxdepth, minweight); if (octname != NULL) fputs(octname, ambfp); fputc('\n', ambfp); @@ -1036,7 +686,6 @@ avstore( /* allocate memory and save aval */ if ((av = newambval()) == NULL) error(SYSTEM, "out of memory in avstore"); *av = *aval; - av->latick = ambclock; av->next = NULL; nambvals++; d = bright(av->val); @@ -1109,178 +758,24 @@ unloadatree( /* unload an ambient value tree */ } -static struct avl { - AMBVAL *p; - unsigned long t; -} *avlist1; /* ambient value list with ticks */ -static AMBVAL **avlist2; /* memory positions for sorting */ -static int i_avlist; /* index for lists */ - -static int alatcmp(const void *av1, const void *av2); - static void avfree(AMBVAL *av) { free(av); } + static void -av2list( - AMBVAL *av -) +sortambvals(void) /* resort ambient values */ { -#ifdef DEBUG - if (i_avlist >= nambvals) - error(CONSISTENCY, "too many ambient values in av2list1"); -#endif - avlist1[i_avlist].p = avlist2[i_avlist] = (AMBVAL*)av; - avlist1[i_avlist++].t = av->latick; -} + AMBTREE oldatrunk = atrunk; - -static int -alatcmp( /* compare ambient values for MRA */ - const void *av1, - const void *av2 -) -{ - long lc = ((struct avl *)av2)->t - ((struct avl *)av1)->t; - return(lc<0 ? -1 : lc>0 ? 1 : 0); + atrunk.alist = NULL; + atrunk.kid = NULL; + unloadatree(&oldatrunk, avinsert); } -/* GW NOTE 2002/10/3: - * I used to compare AMBVAL pointers, but found that this was the - * cause of a serious consistency error with gcc, since the optimizer - * uses some dangerous trick in pointer subtraction that - * assumes pointers differ by exact struct size increments. - */ -static int -aposcmp( /* compare ambient value positions */ - const void *avp1, - const void *avp2 -) -{ - long diff = *(char * const *)avp1 - *(char * const *)avp2; - if (diff < 0) - return(-1); - return(diff > 0); -} - - -static int -avlmemi( /* find list position from address */ - AMBVAL *avaddr -) -{ - AMBVAL **avlpp; - - avlpp = (AMBVAL **)bsearch(&avaddr, avlist2, - nambvals, sizeof(AMBVAL *), aposcmp); - if (avlpp == NULL) - error(CONSISTENCY, "address not found in avlmemi"); - return(avlpp - avlist2); -} - - -static void -sortambvals( /* resort ambient values */ - int always -) -{ - AMBTREE oldatrunk; - AMBVAL tav, *tap, *pnext; - int i, j; - /* see if it's time yet */ - if (!always && (ambclock++ < lastsort+sortintvl || - nambvals < SORT_THRESH)) - return; - /* - * The idea here is to minimize memory thrashing - * in VM systems by improving reference locality. - * We do this by periodically sorting our stored ambient - * values in memory in order of most recently to least - * recently accessed. This ordering was chosen so that new - * ambient values (which tend to be less important) go into - * higher memory with the infrequently accessed values. - * Since we expect our values to need sorting less - * frequently as the process continues, we double our - * waiting interval after each call. - * This routine is also called by setambacc() with - * the "always" parameter set to 1 so that the ambient - * tree will be rebuilt with the new accuracy parameter. - */ - if (tracktime) { /* allocate pointer arrays to sort */ - avlist2 = (AMBVAL **)malloc(nambvals*sizeof(AMBVAL *)); - avlist1 = (struct avl *)malloc(nambvals*sizeof(struct avl)); - } else { - avlist2 = NULL; - avlist1 = NULL; - } - if (avlist1 == NULL) { /* no time tracking -- rebuild tree? */ - if (avlist2 != NULL) - free(avlist2); - if (always) { /* rebuild without sorting */ - oldatrunk = atrunk; - atrunk.alist = NULL; - atrunk.kid = NULL; - unloadatree(&oldatrunk, avinsert); - } - } else { /* sort memory by last access time */ - /* - * Sorting memory is tricky because it isn't contiguous. - * We have to sort an array of pointers by MRA and also - * by memory position. We then copy values in "loops" - * to minimize memory hits. Nevertheless, we will visit - * everyone at least twice, and this is an expensive process - * when we're thrashing, which is when we need to do it. - */ -#ifdef DEBUG - sprintf(errmsg, "sorting %u ambient values at ambclock=%lu...", - nambvals, ambclock); - eputs(errmsg); -#endif - i_avlist = 0; - unloadatree(&atrunk, av2list); /* empty current tree */ -#ifdef DEBUG - if (i_avlist < nambvals) - error(CONSISTENCY, "missing ambient values in sortambvals"); -#endif - qsort(avlist1, nambvals, sizeof(struct avl), alatcmp); - qsort(avlist2, nambvals, sizeof(AMBVAL *), aposcmp); - for (i = 0; i < nambvals; i++) { - if (avlist1[i].p == NULL) - continue; - tap = avlist2[i]; - tav = *tap; - for (j = i; (pnext = avlist1[j].p) != tap; - j = avlmemi(pnext)) { - *(avlist2[j]) = *pnext; - avinsert(avlist2[j]); - avlist1[j].p = NULL; - } - *(avlist2[j]) = tav; - avinsert(avlist2[j]); - avlist1[j].p = NULL; - } - free(avlist1); - free(avlist2); - /* compute new sort interval */ - sortintvl = ambclock - lastsort; - if (sortintvl >= MAX_SORT_INTVL/2) - sortintvl = MAX_SORT_INTVL; - else - sortintvl <<= 1; /* wait twice as long next */ -#ifdef DEBUG - eputs("done\n"); -#endif - } - if (ambclock >= MAXACLOCK) - ambclock = MAXACLOCK/2; - lastsort = ambclock; -} - - #ifdef F_SETLKW static void @@ -1292,9 +787,14 @@ aflock( /* lock/unlock ambient file */ if (typ == fls.l_type) /* already called? */ return; + fls.l_type = typ; - if (fcntl(fileno(ambfp), F_SETLKW, &fls) < 0) - error(SYSTEM, "cannot (un)lock ambient file"); + do + if (fcntl(fileno(ambfp), F_SETLKW, &fls) != -1) + return; + while (errno == EINTR); + + error(SYSTEM, "cannot (un)lock ambient file"); } @@ -1313,10 +813,10 @@ ambsync(void) /* synchronize ambient file */ if ((flen = lseek(fileno(ambfp), (off_t)0, SEEK_END)) < 0) goto seekerr; if ((n = flen - lastpos) > 0) { /* file has grown */ - if (ambinp == NULL) { /* use duplicate filedes */ - ambinp = fdopen(dup(fileno(ambfp)), "r"); + if (ambinp == NULL) { /* get new file pointer */ + ambinp = fopen(ambfile, "rb"); if (ambinp == NULL) - error(SYSTEM, "fdopen failed in ambsync"); + error(SYSTEM, "fopen failed in ambsync"); } if (fseek(ambinp, lastpos, SEEK_SET) < 0) goto seekerr; @@ -1331,24 +831,18 @@ ambsync(void) /* synchronize ambient file */ avstore(&avs); n -= AMBVALSIZ; } - lastpos = flen - n; - /*** seek always as safety measure - if (n) ***/ /* alignment */ - if (lseek(fileno(ambfp), (off_t)lastpos, SEEK_SET) < 0) - goto seekerr; + lastpos = flen - n; /* check alignment */ + if (n && lseek(fileno(ambfp), (off_t)lastpos, SEEK_SET) < 0) + goto seekerr; } n = fflush(ambfp); /* calls write() at last */ - if (n != EOF) - lastpos += (long)nunflshed*AMBVALSIZ; - else if ((lastpos = lseek(fileno(ambfp), (off_t)0, SEEK_CUR)) < 0) - goto seekerr; - + lastpos += (long)nunflshed*AMBVALSIZ; aflock(F_UNLCK); /* release file */ nunflshed = 0; return(n); seekerr: error(SYSTEM, "seek failed in ambsync"); - return -1; /* pro forma return */ + return(EOF); /* pro forma return */ } #else /* ! F_SETLKW */