/** * Clip a ray against a rectangular parallelepiped (RPP) that has faces * parallel to the coordinate planes (a clipping RPP). The RPP is * defined by a minimum point and a maximum point. * * Returns - * 0 if ray does not hit RPP, * !0 if ray hits RPP. * * Implicit Return - * if !0 was returned, "a" and "b" have been clipped to the RPP. */ int vclip(vect_t a, vect_t b, fastf_t *min, fastf_t *max) { static vect_t diff; static double sv; static double st; static double mindist, maxdist; fastf_t *pt = &a[0]; fastf_t *dir = &diff[0]; int i; mindist = -CLIP_DISTANCE; maxdist = CLIP_DISTANCE; VSUB2(diff, b, a); for (i=0; i < 3; i++, pt++, dir++, max++, min++) { if (*dir < -EPSILON) { if ((sv = (*min - *pt) / *dir) < 0.0) return 0; /* MISS */ if (maxdist > sv) maxdist = sv; if (mindist < (st = (*max - *pt) / *dir)) mindist = st; } else if (*dir > EPSILON) { if ((st = (*max - *pt) / *dir) < 0.0) return 0; /* MISS */ if (maxdist > st) maxdist = st; if (mindist < ((sv = (*min - *pt) / *dir))) mindist = sv; } else { /* * If direction component along this axis is NEAR 0, * (i.e., this ray is aligned with this axis), merely * check against the boundaries. */ if ((*min > *pt) || (*max < *pt)) return 0; /* MISS */; } } if (mindist >= maxdist) return 0; /* MISS */ if (mindist > 1 || maxdist < 0) return 0; /* MISS */ if (mindist <= 0 && maxdist >= 1) return 1; /* HIT, no clipping needed */ /* Don't grow one end of a contained segment */ if (mindist < 0) mindist = 0; if (maxdist > 1) maxdist = 1; /* Compute actual intercept points */ VJOIN1(b, a, maxdist, diff); /* b must go first */ VJOIN1(a, a, mindist, diff); return 1; /* HIT */ }
inline void rt_metaball_norm_internal(vect_t *n, point_t *p, struct rt_metaball_internal *mb) { struct wdb_metaballpt *mbpt; vect_t v; fastf_t a; VSETALL(*n, 0.0); switch (mb->method) { case METABALL_METABALL: bu_log("Sorry, strict metaballs are not yet implemented\n"); break; case METABALL_ISOPOTENTIAL: for (BU_LIST_FOR(mbpt, wdb_metaballpt, &mb->metaball_ctrl_head)) { VSUB2(v, *p, mbpt->coord); a = MAGSQ(v); VJOIN1(*n, *n, fabs(mbpt->fldstr)*mbpt->fldstr / (SQ(a)), v); /* f/r^4 */ } break; case METABALL_BLOB: for (BU_LIST_FOR(mbpt, wdb_metaballpt, &mb->metaball_ctrl_head)) { VSUB2(v, *p, mbpt->coord); a = MAGSQ(v); VJOIN1(*n, *n, 2.0*mbpt->sweat/SQ(mbpt->fldstr)*exp(mbpt->sweat*(1-(a/SQ(mbpt->fldstr)))) , v); } break; default: bu_log("unknown metaball method\n"); break; } VUNITIZE(*n); }
HIDDEN int raydiff_hit(struct application *ap, struct partition *PartHeadp, struct seg *UNUSED(segs)) { point_t in_pt, out_pt; struct partition *part; fastf_t part_len = 0.0; struct raydiff_container *state = (struct raydiff_container *)(ap->a_uptr); /*rt_pr_seg(segs);*/ /*rt_pr_partitions(ap->a_rt_i, PartHeadp, "hits");*/ for (part = PartHeadp->pt_forw; part != PartHeadp; part = part->pt_forw) { VJOIN1(in_pt, ap->a_ray.r_pt, part->pt_inhit->hit_dist, ap->a_ray.r_dir); VJOIN1(out_pt, ap->a_ray.r_pt, part->pt_outhit->hit_dist, ap->a_ray.r_dir); part_len = DIST_PT_PT(in_pt, out_pt); if (part_len > state->tol) { state->have_diffs = 1; if (state->left && !bu_strncmp(part->pt_regionp->reg_name+1, state->left_name, strlen(state->left_name))) { RDIFF_ADD_DSEG(state->left, in_pt, out_pt); bu_log("LEFT diff vol (%s) (len: %f): %g %g %g -> %g %g %g\n", part->pt_regionp->reg_name, part_len, V3ARGS(in_pt), V3ARGS(out_pt)); continue; } if (state->right && !bu_strncmp(part->pt_regionp->reg_name+1, state->right_name, strlen(state->right_name))) { RDIFF_ADD_DSEG(state->right, in_pt, out_pt); bu_log("RIGHT diff vol (%s) (len: %f): %g %g %g -> %g %g %g\n", part->pt_regionp->reg_name, part_len, V3ARGS(in_pt), V3ARGS(out_pt)); continue; } /* If we aren't collecting segments, we already have our final answer */ if (!state->left && !state->right) return 0; } } return 0; }
/* add all hit point info to info list */ HIDDEN int add_hit_pnts(struct application *app, struct partition *partH, struct seg *UNUSED(segs)) { struct partition *pp; struct soltab *stp; /*point_t hit_pnt; vect_t hit_normal;*/ struct rt_point_container *c = (struct rt_point_container *)(app->a_uptr); struct npoints *npt; if (c->pnt_cnt > c->capacity-1) { c->capacity *= 4; c->pts = (struct npoints *)bu_realloc((char *)c->pts, c->capacity * sizeof(struct npoints), "enlarge results array"); } RT_CK_APPLICATION(app); /*struct bu_vls *fp = (struct bu_vls *)(app->a_uptr);*/ /* add all hit points */ for (pp = partH->pt_forw; pp != partH; pp = pp->pt_forw) { npt = &(c->pts[c->pnt_cnt]); /* add "in" hit point info */ stp = pp->pt_inseg->seg_stp; /* hack fix for bad tgc surfaces */ if (bu_strncmp("rec", stp->st_meth->ft_label, 3) == 0 || bu_strncmp("tgc", stp->st_meth->ft_label, 3) == 0) { /* correct invalid surface number */ if (pp->pt_inhit->hit_surfno < 1 || pp->pt_inhit->hit_surfno > 3) { pp->pt_inhit->hit_surfno = 2; } if (pp->pt_outhit->hit_surfno < 1 || pp->pt_outhit->hit_surfno > 3) { pp->pt_outhit->hit_surfno = 2; } } VJOIN1(npt->in.p, app->a_ray.r_pt, pp->pt_inhit->hit_dist, app->a_ray.r_dir); RT_HIT_NORMAL(npt->in.n, pp->pt_inhit, stp, &(app->a_ray), pp->pt_inflip); npt->in.is_set = 1; //bu_vls_printf(fp, "%f %f %f %f %f %f\n", hit_pnt[0], hit_pnt[1], hit_pnt[2], hit_normal[0], hit_normal[1], hit_normal[2]); /* add "out" hit point info (unless half-space) */ stp = pp->pt_inseg->seg_stp; if (bu_strncmp("half", stp->st_meth->ft_label, 4) != 0) { VJOIN1(npt->out.p, app->a_ray.r_pt, pp->pt_outhit->hit_dist, app->a_ray.r_dir); RT_HIT_NORMAL(npt->out.n, pp->pt_outhit, stp, &(app->a_ray), pp->pt_outflip); npt->out.is_set = 1; //bu_vls_printf(fp, "%f %f %f %f %f %f\n", hit_pnt[0], hit_pnt[1], hit_pnt[2], hit_normal[0], hit_normal[1], hit_normal[2]); } c->pnt_cnt++; } return 1; }
void Adjust(void) { fastf_t beta, d, len; struct points *ptr; if ( root == NULL ) return; VMOVE( root->p1, root->p ); VMOVE( root->p2, root->p ); ptr = root->next; if ( ptr == NULL ) return; if ( ptr->next == NULL ) { VMOVE( ptr->p1, ptr->p ); VMOVE( ptr->p2, ptr->p ); return; } while ( ptr->next != NULL ) { if ( !torus && !mitre ) { VMOVE( ptr->p1, ptr->p ); VMOVE( ptr->p2, ptr->p ); } else if ( torus ) { /* beta=.5*( pi-ptr->alpha ); d=( MINR+1.0 )*radius*tan( beta ); */ /* dist from new endpts to p2 */ beta = 0.5 * ptr->alpha; d = (MINR + 1.0) * radius / tan( beta ); VJOIN1( ptr->p1, ptr->p, d, ptr->nprev ); VJOIN1( ptr->p2, ptr->p, d, ptr->nnext ); d = sqrt( d*d + (MINR+1.0)*(MINR+1.)*radius*radius ); VJOIN1( ptr->center, ptr->p, d, ptr->nmitre ); } else if ( mitre ) { len = radius/tan(ptr->alpha/2.0); VJOIN1( ptr->p1, ptr->p, -len, ptr->nprev ); VJOIN1( ptr->p2, ptr->p, -len, ptr->nnext ); } ptr = ptr->next; } VMOVE( ptr->p1, ptr->p ); VMOVE( ptr->p2, ptr->p ); }
int BrlCadInterface::hit(application *ap, struct partition *PartHeadp, seg *segs) { register struct partition *pp; register struct hit *hitp; register struct soltab *stp; struct curvature cur; point_t pt; vect_t inormal; vect_t onormal; double curv; int N = 0; //for (pp=PartHeadp->pt_forw; pp != PartHeadp; pp = pp->pt_forw) { pp = PartHeadp->pt_forw; { ++N; hitp = pp->pt_inhit; stp = pp->pt_inseg->seg_stp; VJOIN1(pt, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); RT_HIT_NORMAL(inormal, hitp, stp, &(ap->a_ray), pp->pt_inflip); RT_CURVATURE(&cur, hitp, pp->pt_inflip, stp); curv = max(fabs(cur.crv_c1), fabs(cur.crv_c2)); m_InRadius = 1.0/max(1e-10, curv); m_XIn[0] = pt[0]; m_XIn[1] = pt[1]; m_XIn[2] = pt[2]; m_InNormal[0] = inormal[0]; m_InNormal[1] = inormal[1]; m_InNormal[2] = inormal[2]; hitp = pp->pt_outhit; stp = pp->pt_outseg->seg_stp; VJOIN1(pt, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); RT_HIT_NORMAL( onormal, hitp, stp, &(ap->a_ray), pp->pt_outflip ); RT_CURVATURE(&cur, hitp, pp->pt_inflip, stp); curv = max(fabs(cur.crv_c1), fabs(cur.crv_c2)); m_OutRadius = 1.0/max(1e-10, curv); m_XOut[0] = pt[0]; m_XOut[1] = pt[1]; m_XOut[2] = pt[2]; m_OutNormal[0] = onormal[0]; m_OutNormal[1] = onormal[1]; m_OutNormal[2] = onormal[2]; } m_Hit = true; return(0); }
/** * R E C _ N O R M * * Given ONE ray distance, return the normal and entry/exit point. * hit_surfno is a flag indicating if normal needs to be computed or not. */ void rt_rec_norm(struct hit *hitp, struct soltab *stp, struct xray *rp) { struct rec_specific *rec = (struct rec_specific *)stp->st_specific; VJOIN1(hitp->hit_point, rp->r_pt, hitp->hit_dist, rp->r_dir); switch (hitp->hit_surfno) { case REC_NORM_BODY: /* compute it */ hitp->hit_vpriv[Z] = 0.0; MAT4X3VEC(hitp->hit_normal, rec->rec_invRoS, hitp->hit_vpriv); VUNITIZE(hitp->hit_normal); break; case REC_NORM_TOP: VMOVE(hitp->hit_normal, rec->rec_Hunit); break; case REC_NORM_BOT: VREVERSE(hitp->hit_normal, rec->rec_Hunit); break; default: bu_log("rt_rec_norm: surfno=%d bad\n", hitp->hit_surfno); break; } }
static int hit_headon(struct application *ap, struct partition *PartHeadp) { register char diff_solid; vect_t diff; register fastf_t len; if (PartHeadp->pt_forw->pt_forw != PartHeadp) Tcl_AppendResult(interp, "hit_headon: multiple partitions\n", (char *)NULL); VJOIN1(PartHeadp->pt_forw->pt_inhit->hit_point, ap->a_ray.r_pt, PartHeadp->pt_forw->pt_inhit->hit_dist, ap->a_ray.r_dir); VSUB2(diff, PartHeadp->pt_forw->pt_inhit->hit_point, aim_point); diff_solid = (FIRST_SOLID(sp) != PartHeadp->pt_forw->pt_inseg->seg_stp->st_dp); len = MAGNITUDE(diff); if ( NEAR_ZERO(len, epsilon) || ( diff_solid && VDOT(diff, ap->a_ray.r_dir) > 0 ) ) return(1); else return(0); }
static void plot_ray_img(struct application *ap, const struct partition *pp, double dist, struct bbd_img *bi) { static int plot_num; FILE *pfd; char name[256]; point_t pt; sprintf(name, "bbd_%d.plot3", plot_num++); bu_log("plotting %s\n", name); if ((pfd = fopen(name, "wb")) == (FILE *)NULL) { bu_bomb("can't open plot3 file\n"); } /* red line from ray origin to hit point */ VJOIN1(pp->pt_inhit->hit_point, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir); if (VNEAR_EQUAL(ap->a_ray.r_pt, pp->pt_inhit->hit_point, 0.125)) { /* start and hit point identical, make special allowance */ vect_t vtmp; pl_color(pfd, 255, 0, 128); VREVERSE(vtmp, ap->a_ray.r_dir); VJOIN1(vtmp, ap->a_ray.r_pt, 5.0, vtmp); pdv_3line(pfd, vtmp, pp->pt_inhit->hit_point); } else { pl_color(pfd, 255, 0, 0); pdv_3line(pfd, ap->a_ray.r_pt, pp->pt_inhit->hit_point); } /* yellow line from hit point to plane point */ VJOIN1(pt, ap->a_ray.r_pt, dist, ap->a_ray.r_dir); /* point on plane */ pl_color(pfd, 255, 255, 0); pdv_3line(pfd, pp->pt_inhit->hit_point, pt); /* green line from image origin to plane point */ pl_color(pfd, 0, 255, 0); pdv_3line(pfd, pt, bi->img_origin); fclose(pfd); }
/** * Draw a ray */ void pdv_3ray(FILE *fp, const fastf_t *pt, const fastf_t *dir, double t) { point_t tip; VJOIN1( tip, pt, t, dir ); pdv_3move( fp, pt ); pdv_3cont( fp, tip ); }
/** * Given ONE ray distance, return the normal and entry/exit point. */ void rt_xxx_norm(struct hit *hitp, struct soltab *stp, struct xray *rp) { struct xxx_specific *xxx; if (!stp) return; xxx = (struct xxx_specific *)stp->st_specific; if (!xxx) return; RT_CK_SOLTAB(stp); VJOIN1(hitp->hit_point, rp->r_pt, hitp->hit_dist, rp->r_dir); }
int miss(register struct application *ap) { bu_log("missed\n"); if ( R_DEBUG&RDEBUG_RAYPLOT ) { vect_t out; VJOIN1( out, ap->a_ray.r_pt, 10000, ap->a_ray.r_dir ); /* to imply direction */ pl_color( plotfp, 190, 0, 0 ); pdv_3line( plotfp, ap->a_ray.r_pt, out ); } return(0); }
HIDDEN int raydiff_overlap(struct application *ap, struct partition *pp, struct region *reg1, struct region *reg2, struct partition *UNUSED(hp)) { point_t in_pt, out_pt; fastf_t overlap_len = 0.0; struct raydiff_container *state = (struct raydiff_container *)ap->a_uptr; VJOIN1(in_pt, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir); VJOIN1(out_pt, ap->a_ray.r_pt, pp->pt_outhit->hit_dist, ap->a_ray.r_dir); overlap_len = DIST_PT_PT(in_pt, out_pt); if (overlap_len > state->tol) { RDIFF_ADD_DSEG(state->both, in_pt, out_pt); bu_log("OVERLAP (%s and %s) (len: %f): %g %g %g -> %g %g %g\n", reg1->reg_name, reg2->reg_name, overlap_len, V3ARGS(in_pt), V3ARGS(out_pt)); } return 0; }
/** * Given ONE ray distance, return the normal and entry/exit point. */ void rt_arbn_norm(struct hit *hitp, struct soltab *stp, struct xray *rp) { struct rt_arbn_internal *aip = (struct rt_arbn_internal *)stp->st_specific; size_t h; VJOIN1(hitp->hit_point, rp->r_pt, hitp->hit_dist, rp->r_dir); h = hitp->hit_surfno; if (h > aip->neqn) { bu_log("rt_arbn_norm(): hit_surfno=%zu?\n", h); VSETALL(hitp->hit_normal, 0); return; } VMOVE(hitp->hit_normal, aip->eqn[h]); }
/** * Given ONE ray distance, return the normal and entry/exit point. */ void rt_superell_norm(struct hit *hitp, struct soltab *stp, struct xray *rp) { struct superell_specific *superell = (struct superell_specific *)stp->st_specific; vect_t xlated; fastf_t scale; VJOIN1(hitp->hit_point, rp->r_pt, hitp->hit_dist, rp->r_dir); VSUB2(xlated, hitp->hit_point, superell->superell_V); MAT4X3VEC(hitp->hit_normal, superell->superell_invRSSR, xlated); scale = 1.0 / MAGNITUDE(hitp->hit_normal); VSCALE(hitp->hit_normal, hitp->hit_normal, scale); return; }
void plotRayPoint(struct xray *rayp) { int endpoint[3]; if (plotfp == NULL) return; VJOIN1(endpoint, rayp->r_pt, cellsz, rayp->r_dir); bu_semaphore_acquire(BU_SEM_SYSCALL); pl_color(plotfp, R_BURST, G_BURST, B_BURST); /* draw point */ pl_3point(plotfp, (int) endpoint[X], (int) endpoint[Y], (int) endpoint[Z]); bu_semaphore_release(BU_SEM_SYSCALL); return; }
int rayhit2(struct application *ap, register struct partition *pt, struct seg *UNUSED(segp)) { struct partition *pp = pt->pt_forw; struct hit *hitp = pt->pt_forw->pt_inhit; struct cell *c = (struct cell *)ap->a_uptr; c->c_ishit = 1; c->c_region = pp->pt_regionp; c->c_id = pp->pt_regionp->reg_regionid; VMOVE(c->c_rdir, ap->a_ray.r_dir); VJOIN1(c->c_hit, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); RT_HIT_NORMAL(c->c_normal, hitp, pp->pt_inseg->seg_stp, &(ap->a_ray), pp->pt_inflip); c->c_dist = hitp->hit_dist; return 1; }
/** * Given ONE ray distance, return the normal and entry/exit point. */ void rt_hyp_norm(struct hit *hitp, struct soltab *stp, struct xray *rp) { struct hyp_specific *hyp = (struct hyp_specific *)stp->st_specific; /* normal from basic hyperboloid and transformed normal */ vect_t n, nT; VJOIN1(hitp->hit_point, rp->r_pt, hitp->hit_dist, rp->r_dir); switch (hitp->hit_surfno) { case HYP_NORM_TOP: VMOVE(hitp->hit_normal, hyp->hyp_Hunit); break; case HYP_NORM_BOTTOM: VREVERSE(hitp->hit_normal, hyp->hyp_Hunit); break; case HYP_NORM_BODY: /* normal vector is VUNITIZE(z * dz/dx, z * dz/dy, -z) */ /* z = +- (c/a) * sqrt(x^2/a^2 + y^2/b^2 -1) */ VSET(n, hyp->hyp_rx * hitp->hit_vpriv[X], hyp->hyp_ry * hitp->hit_vpriv[Y], -hyp->hyp_rz * hitp->hit_vpriv[Z]); nT[X] = (hyp->hyp_Aunit[X] * n[X]) + (hyp->hyp_Bunit[X] * n[Y]) + (hyp->hyp_Hunit[X] * n[Z]); nT[Y] = (hyp->hyp_Aunit[Y] * n[X]) + (hyp->hyp_Bunit[Y] * n[Y]) + (hyp->hyp_Hunit[Y] * n[Z]); nT[Z] = (hyp->hyp_Aunit[Z] * n[X]) + (hyp->hyp_Bunit[Z] * n[Y]) + (hyp->hyp_Hunit[Z] * n[Z]); VUNITIZE(nT); VMOVE(hitp->hit_normal, nT); break; default: bu_log("rt_hyp_norm: surfno=%d bad\n", hitp->hit_surfno); break; } }
int raymiss(register struct application *ap) { struct cell *posp; /* store the current cell position */ /* Getting defensive.... just in case. */ if ((size_t)ap->a_x > width) { bu_exit(EXIT_FAILURE, "raymiss: pixels exceed width\n"); } posp = &(cellp[ap->a_x]); /* Find the hit point for the miss. */ VJOIN1(posp->c_hit, ap->a_ray.r_pt, max_dist, ap->a_ray.r_dir); posp->c_dist = max_dist; return 0; }
/* * Given a seg which participates in the partition we are shading evaluate * the transmission on the path */ static double eval_seg(struct application *ap, struct reg_db_internals *dbint, struct seg *seg_p) { double span; point_t pt; struct rt_ell_internal *ell_p = (struct rt_ell_internal *)dbint->ip.idb_ptr; double optical_density = 0.0; double step_dist; double dist; int steps; /* XXX Should map the ray into the coordinate system of the ellipsoid * here, so that all computations are done in an axis-aligned system * with the axes being the gaussian dimensions */ span = seg_p->seg_out.hit_dist - seg_p->seg_in.hit_dist; steps = (int)(span / 100.0 + 0.5); if (steps < 2) steps = 2; step_dist = span / (double)steps; if (rdebug&RDEBUG_SHADE) { bu_log("Evaluating Segment:\n"); bu_log("dist_in:%g dist_out:%g span:%g step_dist:%g steps:%d\n", seg_p->seg_in.hit_dist, seg_p->seg_out.hit_dist, span, step_dist, steps); } for (dist=seg_p->seg_in.hit_dist; dist < seg_p->seg_out.hit_dist; dist += step_dist) { VJOIN1(pt, ap->a_ray.r_pt, dist, ap->a_ray.r_dir); optical_density += gauss_eval(pt, ell_p->v, dbint->one_sigma); } return optical_density; }
int rayhit(struct application *ap, register struct partition *PartHeadp, struct seg *UNUSED(segp)) { register struct partition *pp = PartHeadp->pt_forw; struct cell *posp; /* stores current cell position */ if ( pp == PartHeadp ) return 0; /* nothing was actually hit?? */ /* Getting defensive.... just in case. */ if ((size_t)ap->a_x > width) { bu_exit(EXIT_FAILURE, "rayhit: pixels exceed width\n"); } posp = &(cellp[ap->a_x]); /* Calculate the hit distance and the direction vector. This is done * by VJOIN1(hitp->hit_point, rp->r_pt, hitp->hit_dist, rp->r_dir). */ VJOIN1(pp->pt_inhit->hit_point, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir); /* Now store the distance and the direction vector as appropriate. * Output the ray data: screen plane (pixel) coordinates * for x and y positions of a ray, region_id, and hit_distance. * The x and y positions are represented by ap->a_x and ap->a_y. * * Assume all rays are parallel. */ posp->c_dist = pp->pt_inhit->hit_dist; VMOVE(posp->c_hit, pp->pt_inhit->hit_point); return 0; }
/* * Find which section a given X value is in, and indicate what * direction the tube is headed in then. */ void xfinddir(fastf_t *dir, double x, fastf_t *loc) { size_t i; fastf_t ratio; for (i=0; i<nsamples-1; i++) { if (x < sample[i][X]) break; if (x >= sample[i+1][X]) continue; goto out; } fprintf(stderr, "xfinddir: x=%g is past last segment, using final direction\n", x); i = nsamples-2; out: VSUB2(dir, sample[i+1], sample[i]); ratio = (x-sample[i][X]) / (sample[i+1][X]-sample[i][X]); VJOIN1(loc, sample[i], ratio, dir); VUNITIZE(dir); return; }
/* * R A Y H I T * * Rayhit() is called by rt_shootray() when the ray hits one or more objects. * A per-shotline header record is written, followed by information about * each object hit. * * Note that the GIFT-3 format uses a different convention for the "zero" * distance along the ray. RT has zero at the ray origin (emanation plain), * while GIFT has zero at the screen plain translated so that it contains * the model origin. This difference is compensated for by adding the * 'dcorrection' distance correction factor. * * Also note that the GIFT-3 format requires information about the start * point of the ray in two formats. First, the h, v coordinates of the * grid cell CENTERS (in screen space coordinates) are needed. * Second, the ACTUAL h, v coordinates fired from are needed. * * An optional rtg3.pl UnixPlot file is written, permitting a * color vector display of ray-model intersections. */ int rayhit(struct application *ap, register struct partition *PartHeadp, struct seg *segp) { register struct partition *pp = PartHeadp->pt_forw; int comp_count; /* component count */ fastf_t dfirst, dlast; /* ray distances */ static fastf_t dcorrection = 0; /* RT to GIFT dist corr */ int card_count; /* # comp. on this card */ const char *fmt; /* printf() format string */ struct bu_vls str; char buf[128]; /* temp. sprintf() buffer */ point_t hv; /* GIFT h, v coords, in inches */ point_t hvcen; int prev_id=-1; point_t first_hit; int first; if ( pp == PartHeadp ) return(0); /* nothing was actually hit?? */ if ( ap->a_rt_i->rti_save_overlaps ) rt_rebuild_overlaps( PartHeadp, ap, 1 ); part_compact(ap, PartHeadp, TOL); /* count components in partitions */ comp_count = 0; for ( pp=PartHeadp->pt_forw; pp!=PartHeadp; pp=pp->pt_forw ) { if ( pp->pt_regionp->reg_regionid > 0 ) { prev_id = pp->pt_regionp->reg_regionid; comp_count++; } else if ( prev_id <= 0 ) { /* normally air would be output along with a solid partition, but this will require a '111' partition */ prev_id = pp->pt_regionp->reg_regionid; comp_count++; } else prev_id = pp->pt_regionp->reg_regionid; } pp = PartHeadp->pt_back; if ( pp!=PartHeadp && pp->pt_regionp->reg_regionid <= 0 ) comp_count++; /* a trailing '111' ident */ if ( comp_count == 0 ) return( 0 ); /* Set up variable length string, to buffer this shotline in. * Note that there is one component per card, and that each card * (line) is 80 characters long. Hence the parameters given to * rt-vls-extend(). */ bu_vls_init( &str ); bu_vls_extend( &str, 80 * (comp_count+1) ); /* * Find the H, V coordinates of the grid cell center. * RT uses the lower left corner of each cell. */ { point_t center; fastf_t dx; fastf_t dy; dx = ap->a_x + 0.5; dy = ap->a_y + 0.5; VJOIN2( center, viewbase_model, dx, dx_model, dy, dy_model ); MAT4X3PNT( hvcen, model2hv, center ); } /* * Find exact h, v coordinates of actual ray start by * projecting start point into GIFT h, v coordinates. */ MAT4X3PNT( hv, model2hv, ap->a_ray.r_pt ); /* * In RT, rays are launched from the plane of the screen, * and ray distances are relative to the start point. * In GIFT-3 output files, ray distances are relative to * the (H, V) plane translated so that it contains the origin. * A distance correction is required to convert between the two. * Since this really should be computed only once, not every time, * the trip_count flag was added. */ { static int trip_count; vect_t tmp; vect_t viewZdir; if ( trip_count == 0) { VSET( tmp, 0, 0, -1 ); /* viewing direction */ MAT4X3VEC( viewZdir, view2model, tmp ); VUNITIZE( viewZdir ); /* dcorrection will typically be negative */ dcorrection = VDOT( ap->a_ray.r_pt, viewZdir ); trip_count = 1; } } /* This code is for diagnostics. * bu_log("dcorrection=%g\n", dcorrection); */ /* dfirst and dlast have been made negative to account for GIFT looking * in the opposite direction of RT. */ dfirst = -(PartHeadp->pt_forw->pt_inhit->hit_dist + dcorrection); dlast = -(PartHeadp->pt_back->pt_outhit->hit_dist + dcorrection); #if 0 /* This code is to note any occurances of negative distances. */ if ( PartHeadp->pt_forw->pt_inhit->hit_dist < 0) { bu_log("ERROR: dfirst=%g at partition x%x\n", dfirst, PartHeadp->pt_forw ); bu_log("\tdcorrection = %f\n", dcorrection ); bu_log("\tray start point is ( %f %f %f ) in direction ( %f %f %f )\n", V3ARGS( ap->a_ray.r_pt ), V3ARGS( ap->a_ray.r_dir ) ); VJOIN1( PartHeadp->pt_forw->pt_inhit->hit_point, ap->a_ray.r_pt, PartHeadp->pt_forw->pt_inhit->hit_dist, ap->a_ray.r_dir ); VJOIN1( PartHeadp->pt_back->pt_outhit->hit_point, ap->a_ray.r_pt, PartHeadp->pt_forw->pt_outhit->hit_dist, ap->a_ray.r_dir ); rt_pr_partitions(ap->a_rt_i, PartHeadp, "Defective partion:"); } /* End of bug trap. */ #endif /* * Output the ray header. The GIFT statements that * would have generated this are: * 410 write(1, 411) hcen, vcen, h, v, ncomp, dfirst, dlast, a, e * 411 format(2f7.1, 2f9.3, i3, 2f8.2,' A', f6.1,' E', f6.1) */ #define SHOT_FMT "%7.1f%7.1f%9.3f%9.3f%3d%8.2f%8.2f A%6.1f E%6.1f" if ( rt_perspective > 0 ) { bn_ae_vec( &azimuth, &elevation, ap->a_ray.r_dir ); } bu_vls_printf( &str, SHOT_FMT, hvcen[0], hvcen[1], hv[0], hv[1], comp_count, dfirst * MM2IN, dlast * MM2IN, azimuth, elevation ); /* * As an aid to debugging, take advantage of the fact that * there are more than 80 columns on UNIX "cards", and * add debugging information to the end of the line to * allow this shotline to be reproduced offline. * -b gives the shotline x, y coordinates when re-running RTG3, * -p and -d are used with RTSHOT * The easy way to activate this is with the harmless -!1 option * when running RTG3. */ if ( R_DEBUG || bu_debug || RT_G_DEBUG ) { bu_vls_printf( &str, " -b%d,%d -p %26.20e %26.20e %26.20e -d %26.20e %26.20e %26.20e\n", ap->a_x, ap->a_y, V3ARGS(ap->a_ray.r_pt), V3ARGS(ap->a_ray.r_dir) ); } else { bu_vls_putc( &str, '\n' ); } /* loop here to deal with individual components */ card_count = 0; prev_id = -1; first = 1; for ( pp=PartHeadp->pt_forw; pp!=PartHeadp; pp=pp->pt_forw ) { /* * The GIFT statements that would have produced * this output are: * do 632 i=icomp, iend * if (clos(icomp).gt.999.99.or.slos(i).gt.999.9) goto 635 * 632 continue * write(1, 633)(item(i), clos(i), cangi(i), cango(i), * & kspac(i), slos(i), i=icomp, iend) * 633 format(1x, 3(i4, f6.2, 2f5.1, i1, f5.1)) * goto 670 * 635 write(1, 636)(item(i), clos(i), cangi(i), cango(i), * & kspac(i), slos(i), i=icomp, iend) * 636 format(1x, 3(i4, f6.1, 2f5.1, i1, f5.0)) */ fastf_t comp_thickness; /* component line of sight thickness */ fastf_t in_obliq; /* in obliquity angle */ fastf_t out_obliq; /* out obliquity angle */ int region_id; /* solid region's id */ int air_id; /* air id */ fastf_t dot_prod; /* dot product of normal and ray dir */ fastf_t air_thickness; /* air line of sight thickness */ vect_t normal; /* surface normal */ register struct partition *nextpp = pp->pt_forw; region_id = pp->pt_regionp->reg_regionid; if ( region_id <= 0 && prev_id > 0 ) { /* air region output with previous partition */ prev_id = region_id; continue; } comp_thickness = pp->pt_outhit->hit_dist - pp->pt_inhit->hit_dist; /* The below code is meant to catch components with zero or * negative thicknesses. This is not supposed to be possible, * but the condition has been seen. */ #if 0 if ( comp_thickness <= 0 ) { VJOIN1( pp->pt_inhit->hit_point, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir ); VJOIN1( pp->pt_outhit->hit_point, ap->a_ray.r_pt, pp->pt_outhit->hit_dist, ap->a_ray.r_dir ); bu_log("ERROR: comp_thickness=%g for region id = %d at h=%g, v=%g (x=%d, y=%d), partition at x%x\n", comp_thickness, region_id, hv[0], hv[1], ap->a_x, ap->a_y, pp ); rt_pr_partitions(ap->a_rt_i, PartHeadp, "Defective partion:"); bu_log("Send this output to the BRL-CAD Developers ([email protected])\n"); if ( ! (RT_G_DEBUG & DEBUG_ARB8)) { rt_g.debug |= DEBUG_ARB8; rt_shootray(ap); rt_g.debug &= ~DEBUG_ARB8; } } #endif if ( nextpp == PartHeadp ) { if ( region_id <= 0 ) { /* last partition is air, need a 111 'phantom armor' before AND after */ bu_log( "WARNING: adding 'phantom armor' (id=111) with zero thickness before and after air region %s\n", pp->pt_regionp->reg_name ); region_id = 111; air_id = pp->pt_regionp->reg_aircode; air_thickness = comp_thickness; comp_thickness = 0.0; } else { /* Last partition, no air follows, use code 9 */ air_id = 9; air_thickness = 0.0; } } else if ( region_id <= 0 ) { /* air region, need a 111 'phantom armor' */ bu_log( "WARNING: adding 'phantom armor' (id=111) with zero thickness before air region %s\n", pp->pt_regionp->reg_name ); prev_id = region_id; region_id = 111; air_id = pp->pt_regionp->reg_aircode; air_thickness = comp_thickness; comp_thickness = 0.0; } else if ( nextpp->pt_regionp->reg_regionid <= 0 && nextpp->pt_regionp->reg_aircode != 0 ) { /* Next partition is air region */ air_id = nextpp->pt_regionp->reg_aircode; air_thickness = nextpp->pt_outhit->hit_dist - nextpp->pt_inhit->hit_dist; prev_id = air_id; } else { /* 2 solid regions, maybe with gap */ air_id = 0; air_thickness = nextpp->pt_inhit->hit_dist - pp->pt_outhit->hit_dist; if ( air_thickness < 0.0 ) air_thickness = 0.0; if ( !NEAR_ZERO( air_thickness, 0.1 ) ) { air_id = 1; /* air gap */ if ( R_DEBUG & RDEBUG_HITS ) bu_log("air gap added\n"); } else { air_thickness = 0.0; } prev_id = region_id; } /* * Compute the obliquity angles in degrees, ie, * the "declension" angle down off the normal vector. * RT normals always point outwards; * the "inhit" normal points opposite the ray direction, * the "outhit" normal points along the ray direction. * Hence the one sign change. * XXX this should probably be done with atan2() */ if ( first ) { first = 0; VJOIN1( first_hit, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir ); } out: RT_HIT_NORMAL( normal, pp->pt_inhit, pp->pt_inseg->seg_stp, &(ap->a_ray), pp->pt_inflip ); dot_prod = VDOT( ap->a_ray.r_dir, normal ); if ( dot_prod > 1.0 ) dot_prod = 1.0; if ( dot_prod < -1.0 ) dot_prod = (-1.0); in_obliq = acos( -dot_prod ) * bn_radtodeg; RT_HIT_NORMAL( normal, pp->pt_outhit, pp->pt_outseg->seg_stp, &(ap->a_ray), pp->pt_outflip ); dot_prod = VDOT( ap->a_ray.r_dir, normal ); if ( dot_prod > 1.0 ) dot_prod = 1.0; if ( dot_prod < -1.0 ) dot_prod = (-1.0); out_obliq = acos( dot_prod ) * bn_radtodeg; /* Check for exit obliquties greater than 90 degrees. */ #if 0 if ( in_obliq > 90 || in_obliq < 0 ) { bu_log("ERROR: in_obliquity=%g\n", in_obliq); rt_pr_partitions(ap->a_rt_i, PartHeadp, "Defective partion:"); } if ( out_obliq > 90 || out_obliq < 0 ) { bu_log("ERROR: out_obliquity=%g\n", out_obliq); VPRINT(" r_dir", ap->a_ray.r_dir); VPRINT("normal", normal); bu_log("dot=%g, acos(dot)=%g\n", VDOT( ap->a_ray.r_dir, normal ), acos( VDOT( ap->a_ray.r_dir, normal ) ) ); /* Print the defective one */ rt_pr_pt( ap->a_rt_i, pp ); /* Print the whole ray's partition list */ rt_pr_partitions(ap->a_rt_i, PartHeadp, "Defective partion:"); } #endif if ( in_obliq > 90.0 ) in_obliq = 90.0; if ( in_obliq < 0.0 ) in_obliq = 0.0; if ( out_obliq > 90.0 ) out_obliq = 90.0; if ( out_obliq < 0.0 ) out_obliq = 0.0; /* * Handle 3-components per card output format, with * a leading space in front of the first component. */ if ( card_count == 0 ) { bu_vls_strcat( &str, " " ); } comp_thickness *= MM2IN; /* Check thickness fields for format overflow */ if ( comp_thickness > 999.99 || air_thickness*MM2IN > 999.9 ) fmt = "%4d%6.1f%5.1f%5.1f%1d%5.0f"; else fmt = "%4d%6.2f%5.1f%5.1f%1d%5.1f"; #ifdef SPRINTF_NOT_PARALLEL bu_semaphore_acquire( BU_SEM_SYSCALL ); #endif snprintf(buf, 128, fmt, region_id, comp_thickness, in_obliq, out_obliq, air_id, air_thickness*MM2IN ); #ifdef SPRINTF_NOT_PARALLEL bu_semaphore_release( BU_SEM_SYSCALL ); #endif bu_vls_strcat( &str, buf ); card_count++; if ( card_count >= 3 ) { bu_vls_strcat( &str, "\n" ); card_count = 0; } /* A color rtg3.pl UnixPlot file of output commands * is generated. This is processed by plot(1) * plotting filters such as pl-fb or pl-sgi. * Portions of a ray passing through air within the * model are represented in blue, while portions * passing through a solid are assigned green. * This will always be done single CPU, * to prevent output garbling. (See view_init). */ if (R_DEBUG & RDEBUG_RAYPLOT) { vect_t inpt; vect_t outpt; VJOIN1(inpt, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir); VJOIN1(outpt, ap->a_ray.r_pt, pp->pt_outhit->hit_dist, ap->a_ray.r_dir); pl_color(plotfp, 0, 255, 0); /* green */ pdv_3line(plotfp, inpt, outpt); if (air_thickness > 0) { vect_t air_end; VJOIN1(air_end, ap->a_ray.r_pt, pp->pt_outhit->hit_dist + air_thickness, ap->a_ray.r_dir); pl_color(plotfp, 0, 0, 255); /* blue */ pdv_3cont(plotfp, air_end); } } if ( nextpp == PartHeadp && air_id != 9 ) { /* need to output a 111 'phantom armor' at end of shotline */ air_id = 9; air_thickness = 0.0; region_id = 111; comp_thickness = 0.0; goto out; } } /* If partway through building the line, add a newline */ if ( card_count > 0 ) { /* * Note that GIFT zero-fills the unused component slots, * but neither COVART II nor COVART III require it, * so just end the line here. */ bu_vls_strcat( &str, "\n" ); } /* Single-thread through file output. * COVART will accept non-sequential ray data provided the * ray header and its associated data are not separated. CAVEAT: * COVART will not accept headers out of sequence. */ bu_semaphore_acquire( BU_SEM_SYSCALL ); fputs( bu_vls_addr( &str ), outfp ); if ( shot_fp ) { fprintf( shot_fp, "%.5f %.5f %.5f %.5f %.5f %.5f %.5f %.5f %ld %.5f %.5f %.5f\n", azimuth, elevation, V3ARGS( ap->a_ray.r_pt ), V3ARGS( ap->a_ray.r_dir ), line_num, V3ARGS( first_hit) ); line_num += 1 + (comp_count / 3 ); if ( comp_count % 3 ) line_num++; } /* End of single-thread region */ bu_semaphore_release( BU_SEM_SYSCALL ); /* Release vls storage */ bu_vls_free( &str ); return(0); }
//******************************************************************* // NAME hit_func // // INPUTS // ap - application structure // PartHeadp - pointer to the beginning of the trace list // // DESCRIPTION // Called when a shotline produces a hit. This user defined // callback copies the in and out hit points into structures and // stores them in the ray's result list. //******************************************************************* static int hit_func(struct application *ap, struct partition *PartHeadp, struct seg *) { /* iterating over partitions, this will keep track of the current * partition we're working on. */ struct partition *pp; /* will serve as a pointer for the entry and exit hitpoints */ struct hit *hitp; /* will serve as a pointer to the solid primitive we hit */ struct soltab *stp; /* output variables */ Interface::ray_results *result = (Interface::ray_results *) ap->a_uptr; Interface::one_hit singleHit; VMOVE(result->origin, ap->a_ray.r_pt); result->x = ap->a_x; result->y = ap->a_y; /* iterate over each partition until we get back to the head. * each partition corresponds to a specific homogeneous region of * material. */ for (pp=PartHeadp->pt_forw; pp != PartHeadp; pp = pp->pt_forw) { /* entry hit point, so we type less */ hitp = pp->pt_inhit; /* primitive we encountered on entry */ stp = pp->pt_inseg->seg_stp; /* clear data struct */ ONE_HIT_INIT(singleHit); /* record IN hit point */ VJOIN1(singleHit.point, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); /* calculate IN normal vector */ /* compute the normal vector at the entry point, flipping the * normal if necessary. */ RT_HIT_NORMAL(singleHit.normal, hitp, stp, &(ap->a_ray), pp->pt_inflip); /* save IN solid and surface numbers */ /* singleHit.solidNumber = stp->st_bit; */ singleHit.solidName = stp->st_dp->d_namep; singleHit.solidSurface = hitp->hit_surfno; /* push hit to the result list */ result->hitData.push_back(singleHit); /* exit point, so we type less */ hitp = pp->pt_outhit; /* primitive we exited from */ stp = pp->pt_outseg->seg_stp; /* clear data struct */ ONE_HIT_INIT(singleHit); /* record OUT hit point */ VJOIN1(singleHit.point, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); /* calculate OUT normal vector */ /* compute the normal vector at the exit point, flipping the * normal if necessary. */ RT_HIT_NORMAL(singleHit.normal, hitp, stp, &(ap->a_ray), pp->pt_outflip); /* save OUT solid and surface numbers */ /* singleHit.solidNumber = stp->st_bit; */ singleHit.solidName = stp->st_dp->d_namep; singleHit.solidSurface = hitp->hit_surfno; /* push hit to the result list */ result->hitData.push_back(singleHit); } return 1; }
/** * rt_shootray() was told to call this on a hit. * * This callback routine utilizes the application structure which * describes the current state of the raytrace. * * This callback routine is provided a circular linked list of * partitions, each one describing one in and out segment of one * region for each region encountered. * * The 'segs' segment list is unused in this example. */ int hit(struct application *ap, struct partition *PartHeadp, struct seg *segs) { /* iterating over partitions, this will keep track of the current * partition we're working on. */ struct partition *pp; /* will serve as a pointer for the entry and exit hitpoints */ struct hit *hitp; /* will serve as a pointer to the solid primitive we hit */ struct soltab *stp; /* will contain surface curvature information at the entry */ struct curvature cur; /* will contain our hit point coordinate */ point_t pt; /* will contain normal vector where ray enters geometry */ vect_t inormal; /* will contain normal vector where ray exits geometry */ vect_t onormal; /* iterate over each partition until we get back to the head. * each partition corresponds to a specific homogeneous region of * material. */ for (pp=PartHeadp->pt_forw; pp != PartHeadp; pp = pp->pt_forw) { /* print the name of the region we hit as well as the name of * the primitives encountered on entry and exit. */ bu_log("\n--- Hit region %s (in %s, out %s)\n", pp->pt_regionp->reg_name, pp->pt_inseg->seg_stp->st_name, pp->pt_outseg->seg_stp->st_name ); /* entry hit point, so we type less */ hitp = pp->pt_inhit; /* construct the actual (entry) hit-point from the ray and the * distance to the intersection point (i.e., the 't' value). */ VJOIN1(pt, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); /* primitive we encountered on entry */ stp = pp->pt_inseg->seg_stp; /* compute the normal vector at the entry point, flipping the * normal if necessary. */ RT_HIT_NORMAL(inormal, hitp, stp, &(ap->a_ray), pp->pt_inflip); /* print the entry hit point info */ rt_pr_hit(" In", hitp); VPRINT( " Ipoint", pt); VPRINT( " Inormal", inormal); /* This next macro fills in the curvature information which * consists on a principle direction vector, and the inverse * radii of curvature along that direction and perpendicular * to it. Positive curvature bends toward the outward * pointing normal. */ RT_CURVATURE(&cur, hitp, pp->pt_inflip, stp); /* print the entry curvature information */ VPRINT("PDir", cur.crv_pdir); bu_log(" c1=%g\n", cur.crv_c1); bu_log(" c2=%g\n", cur.crv_c2); /* exit point, so we type less */ hitp = pp->pt_outhit; /* construct the actual (exit) hit-point from the ray and the * distance to the intersection point (i.e., the 't' value). */ VJOIN1(pt, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); /* primitive we exited from */ stp = pp->pt_outseg->seg_stp; /* compute the normal vector at the exit point, flipping the * normal if necessary. */ RT_HIT_NORMAL(onormal, hitp, stp, &(ap->a_ray), pp->pt_outflip); /* print the exit hit point info */ rt_pr_hit(" Out", hitp); VPRINT( " Opoint", pt); VPRINT( " Onormal", onormal); } /* A more complicated application would probably fill in a new * local application structure and describe, for example, a * reflected or refracted ray, and then call rt_shootray() for * those rays. */ /* Hit routine callbacks generally return 1 on hit or 0 on miss. * This value is returned by rt_shootray(). */ return 1; }
/* Color pixel based on the energy of a point light source (Eps) plus some diffuse illumination (Epd) reflected from the point <x, y> : E = Epd + Eps (1) The energy reflected from diffuse illumination is the product of the reflectance coefficient at point P (Rp) and the diffuse illumination (Id) : Epd = Rp * Id (2) The energy reflected from the point light source is calculated by the sum of the diffuse reflectance (Rd) and the specular reflectance (Rs), multiplied by the intensity of the light source (Ips) : Eps = (Rd + Rs) * Ips (3) The diffuse reflectance is calculated by the product of the reflectance coefficient (Rp) and the cosine of the angle of incidence (I) : Rd = Rp * cos(I) (4) The specular reflectance is calculated by the product of the specular reflectance coefficient and (the cosine of the angle (S) raised to the nth power) : Rs = W(I) * cos(S)**n (5) Where, I is the angle of incidence. S is the angle between the reflected ray and the observer. W returns the specular reflection coefficient as a function of the angle of incidence. n (roughly 1 to 10) represents the shininess of the surface. * This is the heart of the lighting model which is based on a model developed by Bui-Tuong Phong, [see Wm M. Newman and R. F. Sproull, "Principles of Interactive Computer Graphics", McGraw-Hill, 1979] Er = Ra(m)*cos(Ia) + Rd(m)*cos(I1) + W(I1, m)*cos(s)^^n where, Er is the energy reflected in the observer's direction. Ra is the diffuse reflectance coefficient at the point of intersection due to ambient lighting. Ia is the angle of incidence associated with the ambient light source (angle between ray direction (negated) and surface normal). Rd is the diffuse reflectance coefficient at the point of intersection due to primary lighting. I1 is the angle of incidence associated with the primary light source (angle between light source direction and surface normal). m is the material identification code. W is the specular reflectance coefficient, a function of the angle of incidence, range 0.0 to 1.0, for the material. s is the angle between the reflected ray and the observer. ` n 'Shininess' of the material, range 1 to 10. */ HIDDEN int phong_render(register struct application *ap, const struct partition *pp, struct shadework *swp, void *dp) { struct light_specific *lp; #ifndef RT_MULTISPECTRAL fastf_t *intensity; fastf_t dist; point_t pt; vect_t color; #endif fastf_t *to_light; fastf_t cosine; fastf_t refl; int i; vect_t reflected; vect_t work; #ifdef RT_MULTISPECTRAL struct bn_tabdata *ms_matcolor = BN_TABDATA_NULL; #else point_t matcolor; /* Material color */ #endif struct phong_specific *ps = (struct phong_specific *)dp; if (!ps || ps->magic != PL_MAGIC) bu_bomb("phong_render: bad magic\n"); if (pp == NULL) bu_bomb("phong_render: bad partition\n"); if (rdebug&RDEBUG_SHADE) bu_struct_print("phong_render", phong_parse, (char *)ps); swp->sw_transmit = ps->transmit; swp->sw_reflect = ps->reflect; swp->sw_refrac_index = ps->refrac_index; swp->sw_extinction = ps->extinction; #if SW_SET_TRANSMIT if (swp->sw_phong_set_vector & SW_SET_TRANSMIT) swp->sw_transmit = swp->sw_phong_transmit; if (swp->sw_phong_set_vector & SW_SET_REFLECT) swp->sw_reflect = swp->sw_phong_reflect; if (swp->sw_phong_set_vector & SW_SET_REFRAC_INDEX) swp->sw_refrac_index = swp->sw_phong_ri; if (swp->sw_phong_set_vector & SW_SET_EXTINCTION) swp->sw_extinction = swp->sw_phong_extinction; #endif /* SW_SET_TRANSMIT */ if (swp->sw_xmitonly) { if (swp->sw_xmitonly > 1) return 1; /* done -- wanted parameters only */ if (swp->sw_reflect > 0 || swp->sw_transmit > 0) { if (rdebug&RDEBUG_SHADE) bu_log("calling rr_render from phong, sw_xmitonly\n"); (void)rr_render(ap, pp, swp); } return 1; /* done */ } #ifdef RT_MULTISPECTRAL ms_matcolor = bn_tabdata_dup(swp->msw_color); #else VMOVE(matcolor, swp->sw_color); #endif /* Photon Mapping */ #ifndef RT_MULTISPECTRAL color[0]= swp->sw_color[0]; color[1]= swp->sw_color[1]; color[2]= swp->sw_color[2]; #endif #ifndef RT_MULTISPECTRAL if (!PM_Visualize) #endif { /* Diffuse reflectance from "Ambient" light source (at eye) */ if ((cosine = -VDOT(swp->sw_hit.hit_normal, ap->a_ray.r_dir)) > 0.0) { if (cosine > 1.00001) { bu_log("cosAmb=1+%g %s surfno=%d (x%d, y%d, lvl%d)\n", cosine-1, pp->pt_inseg->seg_stp->st_dp->d_namep, swp->sw_hit.hit_surfno, ap->a_x, ap->a_y, ap->a_level); VPRINT(" normal", swp->sw_hit.hit_normal); VPRINT(" r_dir ", ap->a_ray.r_dir); cosine = 1; } #if SW_SET_TRANSMIT if (swp->sw_phong_set_vector & SW_SET_AMBIENT) { cosine *= swp->sw_phong_ambient; } else { cosine *= AmbientIntensity; } #else cosine *= AmbientIntensity; #endif #ifdef RT_MULTISPECTRAL bn_tabdata_scale(swp->msw_color, ms_matcolor, cosine); #else VSCALE(swp->sw_color, matcolor, cosine); #endif } else { #ifdef RT_MULTISPECTRAL bn_tabdata_constval(swp->msw_color, 0.0); #else VSETALL(swp->sw_color, 0); #endif } /* Emission. 0..1 is normal range, -1..0 sucks light out, like OpenGL */ #ifdef RT_MULTISPECTRAL { float emission[3]; struct bn_tabdata *ms_emission = BN_TABDATA_NULL; VMOVE(emission, ps->emission); #if SW_SET_TRANSMIT if (swp->sw_phong_set_vector & SW_SET_EMISSION) { VSETALL(emission, swp->sw_phong_emission); } #endif /* XXX Really should get a curve at prep, not expand RGB samples */ BN_GET_TABDATA(ms_emission, spectrum); rt_spect_reflectance_rgb(ms_emission, emission); bn_tabdata_add(swp->msw_color, swp->msw_color, ms_emission); bn_tabdata_free(ms_emission); } #else #if SW_SET_TRANSMIT if (swp->sw_phong_set_vector & SW_SET_EMISSION) { vect_t tmp; VSETALL(tmp, swp->sw_phong_emission); VADD2(swp->sw_color, swp->sw_color, tmp); } else { VADD2(swp->sw_color, swp->sw_color, ps->emission); } #else VADD2(swp->sw_color, swp->sw_color, ps->emission); #endif /* SW_SET_TRANSMIT */ #endif /* With the advent of procedural shaders, the caller can no longer * provide us reliable light visibility information. The hit point * may have been changed by another shader in a stack. There is no * way that anyone else can tell us whether lights are visible. */ light_obs(ap, swp, ps->mfp->mf_inputs); /* Consider effects of each light source */ for (i=ap->a_rt_i->rti_nlights-1; i>=0; i--) { if ((lp = (struct light_specific *)swp->sw_visible[i]) == LIGHT_NULL) continue; if (rdebug & RDEBUG_LIGHT) { bu_log("phong_render light=%s lightfract=%g\n", lp->lt_name, swp->sw_lightfract[i]); } /* Light is not shadowed -- add this contribution */ #ifndef RT_MULTISPECTRAL intensity = swp->sw_intensity+3*i; #endif to_light = swp->sw_tolight+3*i; /* Diffuse reflectance from this light source. */ if ((cosine=VDOT(swp->sw_hit.hit_normal, to_light)) > 0.0) { if (cosine > 1.00001) { bu_log("cosI=1+%g (x%d, y%d, lvl%d)\n", cosine-1, ap->a_x, ap->a_y, ap->a_level); cosine = 1; } /* Get Obj Hit Point For Attenuation */ #ifndef RT_MULTISPECTRAL if (PM_Activated) { VJOIN1(pt, ap->a_ray.r_pt, pp->pt_inhit->hit_dist, ap->a_ray.r_dir); dist= sqrt((pt[0]-lp->lt_pos[0])*(pt[0]-lp->lt_pos[0]) + (pt[1]-lp->lt_pos[1])*(pt[1]-lp->lt_pos[1]) + (pt[2]-lp->lt_pos[2])*(pt[2]-lp->lt_pos[2]))/1000.0; dist= (1.0/(0.1 + 1.0*dist + 0.01*dist*dist)); refl= dist * ps->wgt_diffuse * cosine * swp->sw_lightfract[i] * lp->lt_intensity; /* bu_log("pt: [%.3f][%.3f, %.3f, %.3f]\n", dist, pt[0], pt[1], pt[2]);*/ } else #endif { refl= ps->wgt_diffuse * swp->sw_lightfract[i] * cosine * lp->lt_fraction; } #ifdef RT_MULTISPECTRAL bn_tabdata_incr_mul3_scale(swp->msw_color, lp->lt_spectrum, swp->msw_intensity[i], ms_matcolor, refl); #else VELMUL3(work, matcolor, lp->lt_color, intensity); VJOIN1(swp->sw_color, swp->sw_color, refl, work); #endif } /* Calculate specular reflectance. * Reflected ray = (2 * cos(i) * Normal) - Incident ray. * Cos(s) = Reflected ray DOT Incident ray. */ cosine *= 2; VSCALE(work, swp->sw_hit.hit_normal, cosine); VSUB2(reflected, work, to_light); if ((cosine = -VDOT(reflected, ap->a_ray.r_dir)) > 0) { if (cosine > 1.00001) { bu_log("cosS=1+%g (x%d, y%d, lvl%d)\n", cosine-1, ap->a_x, ap->a_y, ap->a_level); cosine = 1; } refl = ps->wgt_specular * swp->sw_lightfract[i] * lp->lt_fraction * #ifdef PHAST_PHONG /* It is unnecessary to compute the actual * exponential here since phong is just a * gross hack. We approximate re: * Graphics Gems IV "A Fast Alternative to * Phong's Specular Model" Pg 385 */ cosine / (ps->shine - ps->shine*cosine + cosine); #else phg_ipow(cosine, ps->shine); #endif /* PHAST_PHONG */ #ifdef RT_MULTISPECTRAL bn_tabdata_incr_mul2_scale(swp->msw_color, lp->lt_spectrum, swp->msw_intensity[i], refl); #else VELMUL(work, lp->lt_color, intensity); VJOIN1(swp->sw_color, swp->sw_color, refl, work); #endif } } #ifndef RT_MULTISPECTRAL if (PM_Activated) { IrradianceEstimate(ap, work, swp->sw_hit.hit_point, swp->sw_hit.hit_normal); VELMUL(work, work, color); VADD2(swp->sw_color, work, swp->sw_color); if (swp->sw_color[0] > 1.0) swp->sw_color[0]= 1.0; if (swp->sw_color[1] > 1.0) swp->sw_color[1]= 1.0; if (swp->sw_color[2] > 1.0) swp->sw_color[2]= 1.0; } } else { if (PM_Activated) { /* IrradianceEstimate(work, swp->sw_hit.hit_point, swp->sw_hit.hit_normal); VELMUL(swp->sw_color, work, color);*/ IrradianceEstimate(ap, swp->sw_color, swp->sw_hit.hit_point, swp->sw_hit.hit_normal); if (swp->sw_color[0] > 1.0) swp->sw_color[0]= 1.0; if (swp->sw_color[1] > 1.0) swp->sw_color[1]= 1.0; if (swp->sw_color[2] > 1.0) swp->sw_color[2]= 1.0; } #endif } if (swp->sw_reflect > 0 || swp->sw_transmit > 0) (void)rr_render(ap, pp, swp); #ifdef RT_MULTISPECTRAL bn_tabdata_free(ms_matcolor); #endif return 1; }
/* * F _ H I D E L I N E */ int f_hideline(ClientData clientData, Tcl_Interp *interp, int argc, char **argv) { FILE *plotfp; char visible; int i, numobjs; char *objname[MAXOBJECTS], title[1]; fastf_t len, u, step; float ratio; vect_t last_move; struct rt_i *rtip; struct resource resource; struct application a; vect_t temp; vect_t last, dir; register struct bn_vlist *vp; CHECK_DBI_NULL; if (argc < 2 || 4 < argc) { struct bu_vls vls; bu_vls_init(&vls); bu_vls_printf(&vls, "help H"); Tcl_Eval(interp, bu_vls_addr(&vls)); bu_vls_free(&vls); return TCL_ERROR; } if ((plotfp = fopen(argv[1], "w")) == NULL) { Tcl_AppendResult(interp, "f_hideline: unable to open \"", argv[1], "\" for writing.\n", (char *)NULL); return TCL_ERROR; } pl_space(plotfp, (int)GED_MIN, (int)GED_MIN, (int)GED_MAX, (int)GED_MAX); /* Build list of objects being viewed */ numobjs = 0; FOR_ALL_SOLIDS(sp) { for (i = 0; i < numobjs; i++) { if ( objname[i] == FIRST_SOLID(sp)->d_namep ) break; } if (i == numobjs) objname[numobjs++] = FIRST_SOLID(sp)->d_namep; } Tcl_AppendResult(interp, "Generating hidden-line drawing of the following regions:\n", (char *)NULL); for (i = 0; i < numobjs; i++) Tcl_AppendResult(interp, "\t", objname[i], "\n", (char *)NULL); /* Initialization for librt */ if ((rtip = rt_dirbuild(dbip->dbi_filename, title, 0)) == RTI_NULL) { Tcl_AppendResult(interp, "f_hideline: unable to open model file \"", dbip->dbi_filename, "\"\n", (char *)NULL); return TCL_ERROR; } a.a_hit = hit_headon; a.a_miss = hit_tangent; a.a_overlap = hit_overlap; a.a_rt_i = rtip; a.a_resource = &resource; a.a_level = 0; a.a_onehit = 1; a.a_diverge = 0; a.a_rbeam = 0; if (argc > 2) { sscanf(argv[2], "%f", &step); step = view_state->vs_Viewscale/step; sscanf(argv[3], "%f", &epsilon); epsilon *= view_state->vs_Viewscale/100; } else { step = view_state->vs_Viewscale/256; epsilon = 0.1*view_state->vs_Viewscale; } for (i = 0; i < numobjs; i++) if (rt_gettree(rtip, objname[i]) == -1) Tcl_AppendResult(interp, "f_hideline: rt_gettree failed on \"", objname[i], "\"\n", (char *)NULL); /* Crawl along the vectors raytracing as we go */ VSET(temp, 0.0, 0.0, -1.0); /* looking at model */ MAT4X3VEC(a.a_ray.r_dir, view_state->vs_view2model, temp); VUNITIZE(a.a_ray.r_dir); FOR_ALL_SOLIDS(sp) { ratio = sp->s_size / VIEWSIZE; /* ignore if small or big */ if (ratio >= dmp->dmr_bound || ratio < 0.001) continue; Tcl_AppendResult(interp, "Primitive\n", (char *)NULL); for ( BU_LIST_FOR( vp, bn_vlist, &(sp->s_vlist) ) ) { register int i; register int nused = vp->nused; register int *cmd = vp->cmd; register point_t *pt = vp->pt; for ( i = 0; i < nused; i++, cmd++, pt++ ) { Tcl_AppendResult(interp, "\tVector\n", (char *)NULL); switch ( *cmd ) { case BN_VLIST_POLY_START: case BN_VLIST_POLY_VERTNORM: break; case BN_VLIST_POLY_MOVE: case BN_VLIST_LINE_MOVE: /* move */ VMOVE(last, *pt); MOVE(last); break; case BN_VLIST_POLY_DRAW: case BN_VLIST_POLY_END: case BN_VLIST_LINE_DRAW: /* setup direction && length */ VSUB2(dir, *pt, last); len = MAGNITUDE(dir); VUNITIZE(dir); visible = FALSE; { struct bu_vls tmp_vls; bu_vls_init(&tmp_vls); bu_vls_printf(&tmp_vls, "\t\tDraw 0 -> %g, step %g\n", len, step); Tcl_AppendResult(interp, bu_vls_addr(&tmp_vls), (char *)NULL); bu_vls_free(&tmp_vls); } for (u = 0; u <= len; u += step) { VJOIN1(aim_point, last, u, dir); MAT4X3PNT(temp, view_state->vs_model2view, aim_point); temp[Z] = 100; /* parallel project */ MAT4X3PNT(a.a_ray.r_pt, view_state->vs_view2model, temp); if (rt_shootray(&a)) { if (!visible) { visible = TRUE; MOVE(aim_point); } } else { if (visible) { visible = FALSE; DRAW(aim_point); } } } if (visible) DRAW(aim_point); VMOVE(last, *pt); /* new last vertex */ } } } } fclose(plotfp); return TCL_OK; }
int f_IR_Model(struct application *ap, Octree *op) { fastf_t octnt_min[3], octnt_max[3]; fastf_t delta = modl_radius / pow_Of_2(ap->a_level); fastf_t point[3] = VINIT_ZERO; /* Intersection point. */ fastf_t norml[3] = VINIT_ZERO; /* Unit normal at point. */ /* Push ray origin along ray direction to intersection point. */ VJOIN1(point, ap->a_ray.r_pt, ap->a_uvec[0], ap->a_ray.r_dir); /* Compute octant RPP. */ octnt_min[X] = op->o_points->c_point[X] - delta; octnt_min[Y] = op->o_points->c_point[Y] - delta; octnt_min[Z] = op->o_points->c_point[Z] - delta; octnt_max[X] = op->o_points->c_point[X] + delta; octnt_max[Y] = op->o_points->c_point[Y] + delta; octnt_max[Z] = op->o_points->c_point[Z] + delta; if (NEAR_EQUAL(point[X], octnt_min[X], epsilon)) /* Intersection point lies on plane whose normal is the negative X-axis. */ { norml[X] = -1.0; norml[Y] = 0.0; norml[Z] = 0.0; } else if (NEAR_EQUAL(point[X], octnt_max[X], epsilon)) /* Intersection point lies on plane whose normal is the positive X-axis. */ { norml[X] = 1.0; norml[Y] = 0.0; norml[Z] = 0.0; } else if (NEAR_EQUAL(point[Y], octnt_min[Y], epsilon)) /* Intersection point lies on plane whose normal is the negative Y-axis. */ { norml[X] = 0.0; norml[Y] = -1.0; norml[Z] = 0.0; } else if (NEAR_EQUAL(point[Y], octnt_max[Y], epsilon)) /* Intersection point lies on plane whose normal is the positive Y-axis. */ { norml[X] = 0.0; norml[Y] = 1.0; norml[Z] = 0.0; } else if (NEAR_EQUAL(point[Z], octnt_min[Z], epsilon)) /* Intersection point lies on plane whose normal is the negative Z-axis. */ { norml[X] = 0.0; norml[Y] = 0.0; norml[Z] = -1.0; } else if (NEAR_EQUAL(point[Z], octnt_max[Z], epsilon)) /* Intersection point lies on plane whose normal is the positive Z-axis. */ { norml[X] = 0.0; norml[Y] = 0.0; norml[Z] = 1.0; } { /* Factor in reflectance from "ambient" light source. */ fastf_t intensity = VDOT(norml, lgts[0].dir); /* Calculate index into false-color table. */ int lgtindex = op->o_temp - ir_min; if (lgtindex > ir_max_index) { bu_log("Temperature (%d) above range of data.\n", op->o_temp); return -1; } if (lgtindex < 0) /* Un-assigned octants get colored grey. */ ap->a_color[0] = ap->a_color[1] = ap->a_color[2] = intensity; else { /* Lookup false-coloring for octant's temperature. */ intensity *= RGB_INVERSE; ap->a_color[0] = (fastf_t) (ir_table[lgtindex][RED]) * intensity; ap->a_color[1] = (fastf_t) (ir_table[lgtindex][GRN]) * intensity; ap->a_color[2] = (fastf_t) (ir_table[lgtindex][BLU]) * intensity; } } return 1; }
/** * R T _ P G _ S H O T * * Function - * Shoot a ray at a polygonal object. * * Returns - * 0 MISS * >0 HIT */ int rt_pg_shot(struct soltab *stp, struct xray *rp, struct application *ap, struct seg *seghead) { struct tri_specific *trip = (struct tri_specific *)stp->st_specific; #define MAXHITS 128 /* # surfaces hit, must be even */ struct hit hits[MAXHITS]; struct hit *hp; size_t nhits; nhits = 0; hp = &hits[0]; /* consider each face */ for (; trip; trip = trip->tri_forw) { fastf_t dn; /* Direction dot Normal */ fastf_t abs_dn; fastf_t k; fastf_t alpha, beta; vect_t wxb; /* vertex - ray_start */ vect_t xp; /* wxb cross ray_dir */ /* * Ray Direction dot N. (N is outward-pointing normal) * wn points inwards, and is not unit length. */ dn = VDOT(trip->tri_wn, rp->r_dir); /* * If ray lies directly along the face, (i.e., dot product * is zero), drop this face. */ abs_dn = dn >= 0.0 ? dn : (-dn); if (abs_dn < SQRT_SMALL_FASTF) continue; VSUB2(wxb, trip->tri_A, rp->r_pt); VCROSS(xp, wxb, rp->r_dir); /* Check for exceeding along the one side */ alpha = VDOT(trip->tri_CA, xp); if (dn < 0.0) alpha = -alpha; if (alpha < 0.0 || alpha > abs_dn) continue; /* Check for exceeding along the other side */ beta = VDOT(trip->tri_BA, xp); if (dn > 0.0) beta = -beta; if (beta < 0.0 || beta > abs_dn) continue; if (alpha+beta > abs_dn) continue; k = VDOT(wxb, trip->tri_wn) / dn; /* For hits other than the first one, might check * to see it this is approx. equal to previous one */ /* If dn < 0, we should be entering the solid. * However, we just assume in/out sorting later will work. * Really should mark and check this! */ VJOIN1(hp->hit_point, rp->r_pt, k, rp->r_dir); /* HIT is within planar face */ hp->hit_magic = RT_HIT_MAGIC; hp->hit_dist = k; VMOVE(hp->hit_normal, trip->tri_N); hp->hit_surfno = trip->tri_surfno; if (++nhits >= MAXHITS) { bu_log("rt_pg_shot(%s): too many hits (%zu)\n", stp->st_name, nhits); break; } hp++; } if (nhits == 0) return 0; /* MISS */ /* Sort hits, Near to Far */ rt_hitsort(hits, nhits); /* Remove duplicate hits. We remove one of a pair of hits when they are 1) close together, and 2) both "entry" or both "exit" occurrences. Two immediate "entry" or two immediate "exit" hits suggest that we hit both of two joined faces, while we want to hit only one. An "entry" followed by an "exit" (or vice versa) suggests that we grazed an edge, and thus we should leave both in the hit list. */ { size_t i, j; for (i=0; i<nhits-1; i++) { fastf_t dist; dist = hits[i].hit_dist - hits[i+1].hit_dist; if (NEAR_ZERO(dist, ap->a_rt_i->rti_tol.dist) && VDOT(hits[i].hit_normal, rp->r_dir) * VDOT(hits[i+1].hit_normal, rp->r_dir) > 0) { for (j=i; j<nhits-1; j++) hits[j] = hits[j+1]; nhits--; i--; } } } if (nhits == 1) nhits = 0; if (nhits&1) { size_t i; static int nerrors = 0; /* message counter */ /* * If this condition exists, it is almost certainly due to * the dn==0 check above. Thus, we will make the last * surface rather thin. * This at least makes the * presence of this solid known. There may be something * better we can do. */ if (nerrors++ < 6) { bu_log("rt_pg_shot(%s): WARNING %zu hits:\n", stp->st_name, nhits); bu_log("\tray start = (%g %g %g) ray dir = (%g %g %g)\n", V3ARGS(rp->r_pt), V3ARGS(rp->r_dir)); for (i=0; i < nhits; i++) { point_t tmp_pt; VJOIN1(tmp_pt, rp->r_pt, hits[i].hit_dist, rp->r_dir); if (VDOT(rp->r_dir, hits[i].hit_normal) < 0.0) bu_log("\tentrance at dist=%f (%g %g %g)\n", hits[i].hit_dist, V3ARGS(tmp_pt)); else bu_log("\texit at dist=%f (%g %g %g)\n", hits[i].hit_dist, V3ARGS(tmp_pt)); } } if (nhits > 2) { fastf_t dot1, dot2; size_t j; /* likely an extra hit, * look for consecutive entrances or exits */ dot2 = 1.0; i = 0; while (i<nhits) { dot1 = dot2; dot2 = VDOT(rp->r_dir, hits[i].hit_normal); if (dot1 > 0.0 && dot2 > 0.0) { /* two consecutive exits, * manufacture an entrance at same distance * as second exit. */ for (j=nhits; j>i; j--) hits[j] = hits[j-1]; /* struct copy */ VREVERSE(hits[i].hit_normal, hits[i].hit_normal); dot2 = VDOT(rp->r_dir, hits[i].hit_normal); nhits++; bu_log("\t\tadding fictitious entry at %f (%s)\n", hits[i].hit_dist, stp->st_name); } else if (dot1 < 0.0 && dot2 < 0.0) { /* two consecutive entrances, * manufacture an exit between them. */ for (j=nhits; j>i; j--) hits[j] = hits[j-1]; /* struct copy */ hits[i] = hits[i-1]; /* struct copy */ VREVERSE(hits[i].hit_normal, hits[i-1].hit_normal); dot2 = VDOT(rp->r_dir, hits[i].hit_normal); nhits++; bu_log("\t\tadding fictitious exit at %f (%s)\n", hits[i].hit_dist, stp->st_name); } i++; } } else { hits[nhits] = hits[nhits-1]; /* struct copy */ VREVERSE(hits[nhits].hit_normal, hits[nhits-1].hit_normal); bu_log("\t\tadding fictitious hit at %f (%s)\n", hits[nhits].hit_dist, stp->st_name); nhits++; } } if (nhits&1) { if (nhits < MAXHITS) { hits[nhits] = hits[nhits-1]; /* struct copy */ VREVERSE(hits[nhits].hit_normal, hits[nhits-1].hit_normal); bu_log("\t\tadding fictitious hit at %f (%s)\n", hits[nhits].hit_dist, stp->st_name); nhits++; } else nhits--; } /* nhits is even, build segments */ { struct seg *segp; size_t i; for (i=0; i < nhits; i += 2) { RT_GET_SEG(segp, ap->a_resource); segp->seg_stp = stp; segp->seg_in = hits[i]; /* struct copy */ segp->seg_out = hits[i+1]; /* struct copy */ BU_LIST_INSERT(&(seghead->l), &(segp->l)); } } return nhits; /* HIT */ }
int handle_main_ray(struct application *ap, register struct partition *PartHeadp, struct seg *segp) { register struct partition *pp; register struct hit *hitp; /* which hit */ struct application a2; struct cell me; struct cell below; struct cell left; struct cell above; struct cell right; double intensity = 1.0; int edge = 0; int cpu; int oc = 1; RGBpixel col; RT_APPLICATION_INIT(&a2); memset(&me, 0, sizeof(struct cell)); memset(&below, 0, sizeof(struct cell)); memset(&left, 0, sizeof(struct cell)); cpu = ap->a_resource->re_cpu; if (PartHeadp == NULL || segp == NULL) { /* The main shotline missed. pack the application struct */ me.c_ishit = 0; me.c_dist = MISS_DIST; me.c_id = MISS_ID; me.c_region = 0; VSETALL(me.c_hit, MISS_DIST); VSETALL(me.c_normal, 0); VMOVE(me.c_rdir, ap->a_ray.r_dir); } else { pp = PartHeadp->pt_forw; hitp = pp->pt_inhit; /* * Stuff the information for this cell. */ me.c_ishit = 1; me.c_id = pp->pt_regionp->reg_regionid; me.c_dist = hitp->hit_dist; me.c_region = pp->pt_regionp; VMOVE(me.c_rdir, ap->a_ray.r_dir); VJOIN1(me.c_hit, ap->a_ray.r_pt, hitp->hit_dist, ap->a_ray.r_dir); RT_HIT_NORMAL(me.c_normal, hitp, pp->pt_inseg->seg_stp, &(ap->a_ray), pp->pt_inflip); } /* * Now, fire a ray for both the cell below and if necessary, the * cell to the left. */ a2.a_hit = rayhit2; a2.a_miss = raymiss2; a2.a_onehit = 1; a2.a_rt_i = ap->a_rt_i; a2.a_resource = ap->a_resource; a2.a_logoverlap = ap->a_logoverlap; VSUB2(a2.a_ray.r_pt, ap->a_ray.r_pt, dy_model); /* below */ VMOVE(a2.a_ray.r_dir, ap->a_ray.r_dir); a2.a_uptr = (void *)&below; rt_shootray(&a2); if (ap->a_x == 0) { /* * For the first pixel in a scanline, we have to shoot to the * left. For each pixel afterword, we save the current cell * info to be used as the left side cell info for the * following pixel */ VSUB2(a2.a_ray.r_pt, ap->a_ray.r_pt, dx_model); /* left */ VMOVE(a2.a_ray.r_dir, ap->a_ray.r_dir); a2.a_uptr = (void *)&left; rt_shootray(&a2); } else { left.c_ishit = saved[cpu]->c_ishit; left.c_id = saved[cpu]->c_id; left.c_dist = saved[cpu]->c_dist; left.c_region = saved[cpu]->c_region; VMOVE(left.c_rdir, saved[cpu]->c_rdir); VMOVE(left.c_hit, saved[cpu]->c_hit); VMOVE(left.c_normal, saved[cpu]->c_normal); } if (both_sides) { VADD2(a2.a_ray.r_pt, ap->a_ray.r_pt, dy_model); /* above */ VMOVE(a2.a_ray.r_dir, ap->a_ray.r_dir); a2.a_uptr = (void *)&above; rt_shootray(&a2); VADD2(a2.a_ray.r_pt, ap->a_ray.r_pt, dx_model); /* right */ VMOVE(a2.a_ray.r_dir, ap->a_ray.r_dir); a2.a_uptr = (void *)&right; rt_shootray(&a2); } /* * Is this pixel an edge? */ if (both_sides) { edge = is_edge(&intensity, ap, &me, &left, &below, &right, &above); } else { edge = is_edge(&intensity, ap, &me, &left, &below, NULL, NULL); } /* * Does this pixel occlude the second geometry? Note that we must * check on edges as well since right side and top edges are * actually misses. */ if (occlusion_mode != OCCLUSION_MODE_NONE) if (me.c_ishit || edge) oc = occludes(ap, &me); /* * Perverse Pixel Painting Paradigm(tm) If a pixel should be * written to the fb, writeable is set. */ if (occlusion_mode == OCCLUSION_MODE_EDGES) writeable[cpu][ap->a_x] = (edge && oc); else if (occlusion_mode == OCCLUSION_MODE_HITS) writeable[cpu][ap->a_x] = ((me.c_ishit || edge) && oc); else if (occlusion_mode == OCCLUSION_MODE_DITHER) { if (edge && oc) writeable[cpu][ap->a_x] = 1; else if (me.c_ishit && oc) { /* * Dither mode. * * For occluding non-edges, only write every other pixel. */ if (oc == 1 && ((ap->a_x + ap->a_y) % 2) == 0) writeable[cpu][ap->a_x] = 1; else if (oc == 2) writeable[cpu][ap->a_x] = 1; else writeable[cpu][ap->a_x] = 0; } else { writeable[cpu][ap->a_x] = 0; } } else { if (edge) writeable[cpu][ap->a_x] = 1; else writeable[cpu][ap->a_x] = 0; } if (edge) { if (both_sides) { choose_color(col, intensity, &me, &left, &below, &right, &above); } else { choose_color(col, intensity, &me, &left, &below, NULL, NULL); } scanline[cpu][ap->a_x*3+RED] = col[RED]; scanline[cpu][ap->a_x*3+GRN] = col[GRN]; scanline[cpu][ap->a_x*3+BLU] = col[BLU]; } else { scanline[cpu][ap->a_x*3+RED] = bgcolor[RED]; scanline[cpu][ap->a_x*3+GRN] = bgcolor[GRN]; scanline[cpu][ap->a_x*3+BLU] = bgcolor[BLU]; } /* * Save the cell info for the next pixel. */ saved[cpu]->c_ishit = me.c_ishit; saved[cpu]->c_id = me.c_id; saved[cpu]->c_dist = me.c_dist; saved[cpu]->c_region = me.c_region; VMOVE(saved[cpu]->c_rdir, me.c_rdir); VMOVE(saved[cpu]->c_hit, me.c_hit); VMOVE(saved[cpu]->c_normal, me.c_normal); return edge; }