Newer
Older
#include <ostream.h>
#include <istream.h>
#include "worldline.h"
// #define CAVE 1 // NewN: use potentially-shared-memory malloc
#include "shmem.h"
#include "specks.h"
#include "partiviewc.h"
#include "kira_parti.h" // declare following as extern "C"
#if !CAVE
# define USE_PLOT 1
#endif
#if USE_PLOT
# include "Plot.H"
#endif
extern struct specklist *kira_get_parti( struct stuff *st, double realtime );
extern int kira_parse_args( struct dyndata *dd, struct stuff *st, int argc, char **argv );
extern int kira_draw( struct dyndata *dd, struct stuff *st, struct specklist *slhead,
Matrix *Tc2w, float radperpix );
extern void kira_HRplot( Fl_Plot *, struct stuff *st, struct dyndata *dd );
extern int kira_HRplot_events( Fl_Plot *plot, int ev );
extern int kira_picked( struct stuff *st, GLuint *hit, struct specklist *sl, int speckno );
extern void kira_maxtrail( struct dyndata *dd, struct stuff *st, int newmax );
enum speckfields {
SPECK_ID = 0, // worldline index (= kira index for single stars,
// unique small negative int for others)
SPECK_MASSLOG = 1, // log10(mass/Msun)
SPECK_NCLUMP = 2, // number of stars in clump
SPECK_TLOG = 3, // log10(Teff)
SPECK_LUM = 4, // L/Lsun?
SPECK_STYPE = 5, // star type: 1ms 2cd 3gs 4bh 5sg 6hb
SPECK_ROOTID = 6, // worldline index of root of clump
SPECK_TREEADDR = 7, // bit-encoded tree address within our clump (0 for isolated stars)
SPECK_RINGSIZE = 8, // size of ring marker
SPECK_MU = 9, // mass ratio; = 0 for leaf nodes
SPECK_SEPVEC = 10, // separation vector[3]
SPECK_NDATAFIELDS = 13
};
static char *fieldnames[] = { // Must match ``enum speckfields'' !!
"id", // worldline index (=? kira index for singles), unique <0 for others
"masslog", // log10(mass)
"nclump", // number of stars in clump, = 1 for singles
"Tlog", // log10( Teff )
"Lum", // L/Lsun?
"stype", // stellar type index
"rootid", // id of root of clump, = our id for singles
"treeaddr", // binary-coded address in clump
"ringsize", // size of ring
"mu", // mass ratio for nonleaf nodes
NULL
};
struct vald {
float min, max, sum;
};
struct trailhead {
int maxtrail, ntrails;
int next; /* ring buffer next-slot-to-use */
struct speck *specks;
};
int ih; // current worldbundle index
real tmin, tmax;
real tcur; // current time
real treq; // requested time
int treearcs; // KIRA_{OFF|ON|CROSS|TICK}
float tickscale; // size of treearc cross mark (frac of sep)
int ringsizer; // KIRA_RINGSEP, KIRA_RINGA
float ringscale; // multiplier for above
float ringmin, ringmax; // range of pixel sizes for ring markers
int tracking, wastracking; // id of particle we're tracking, or zero
Point trackpos; // last known position of tracked particle
float massscale; // scale-factor for masses (e.g. conv to Msun)
int myselseq; // sequence number of ww->sel
int nleafsel;
SelMask *leafsel; // sel bitmasks, indexed by particle id
// selection mapping
SelOp intsrc; // for all particles matching intsrc,
SelOp intdest; // then turn on intdest bit(s).
struct trailhead *trails; // per-star specklist of recent history
int maxtrail;
int maxtrailno;
float trailalpha;
float trailpsize;
int trailonly;
SelOp picksel; // what to do when a star is picked
struct speck *marksp; // current pointer, updated by add_speck
int leafcount; // temp, used in recursion only
int interactsel, unionsel; // ditto
int pickcount; // temp, used in kira_picked
worldstuff( struct stuff *st ) {
this->init( st );
nh = 0;
wh = NULL;
ih = 0;
tmin = 0, tmax = 1;
tcur = treq = 0;
treenodes = KIRA_ON;
treerings = KIRA_OFF;
maxtrail = 50;
trails = NULL;
trails = NULL;
maxtrail = 50;
maxtrailno = 0;
trailonly = 0;
trailalpha = 0.6;
trailpsize = 1.0;
myselseq = 0;
leafsel = NULL;
nleafsel = 0;
parse_selexpr( st, "pick", &picksel, NULL, NULL );
selinit( &trailsel );
selinit( &intsrc );
selinit( &intdest );
extern void kira_add_trail( struct stuff *st, struct worldstuff *, int id, struct speck * );
extern void kira_erase_trail( struct stuff *st, struct worldstuff *, int id );
extern void kira_track_break( struct worldstuff *ww, Point *newtrack );
void kira_invalidate( struct dyndata *dd, struct stuff *st ) {
worldstuff *ww = (worldstuff *)dd->data;
int kira_set( struct dyndata *dd, struct stuff *st, int what, double val )
struct worldstuff *ww = (struct worldstuff *)dd->data;
return 2;
switch(what) {
case KIRA_NODES: ww->treenodes = (int)val; break;
case KIRA_RINGS: ww->treerings = (int)val; break;
case KIRA_TREE: ww->treearcs = (int)val; break;
case KIRA_TICKSCALE: ww->tickscale = val; break;
case KIRA_RINGSIZE: ww->ringsizer = (int)val; break;
case KIRA_RINGSCALE: ww->ringscale = val; break;
case KIRA_RINGMIN: ww->ringmin = val; break;
case KIRA_RINGMAX: ww->ringmax = val; break;
case KIRA_TRACK: if(ww->tracking != (int)val) ww->wastracking = 0;
ww->tracking = (int)val; break;
double kira_get( struct dyndata *dd, struct stuff *st, int what )
struct worldstuff *ww = (struct worldstuff *)dd->data;
case KIRA_NODES: return ww->treenodes;
case KIRA_RINGS: return ww->treerings;
case KIRA_TREE: return ww->treearcs;
case KIRA_TICKSCALE: return ww->tickscale;
case KIRA_RINGSIZE: return ww->ringsizer;
case KIRA_RINGSCALE: return ww->ringscale;
case KIRA_RINGMIN: return ww->ringmin;
case KIRA_RINGMAX: return ww->ringmax;
case KIRA_TRACK: return ww->tracking;
#ifdef NEWSTDIO
return 1;
#else /*old stdio, where filebuf::fd() still existed*/
struct timeval tv;
fd_set fds;
if(s == NULL) return 0;
if(s->rdbuf()->in_avail() || s->eof()) return 1;
if(maxwait < 0) maxwait = 0;
tv.tv_sec = (int)maxwait;
tv.tv_usec = (int) (1000000 * (maxwait - tv.tv_sec));
FD_ZERO(&fds);
FD_SET(s->rdbuf()->fd(), &fds);
return(select(s->rdbuf()->fd() + 1, &fds, NULL, NULL, &tv) > 0);
int kira_help( struct dyndata *dd, struct stuff *st, int verbose )
{
msg(" kira {node|ring|size|scale|span|track|trail} starlab controls; try \"kira ?\"");
return 1;
}
int kira_read_more( struct dyndata *dd, struct stuff *st, int nmax, double tmax, double realdtmax )
struct worldstuff *ww = (struct worldstuff *)dd->data;
if(ww == NULL) return -1;
ifstream *ins = ww->ins;
if(ins == NULL || !ins->rdbuf()->is_open()) return -1;
int paras = 0;
int timelimit = (realdtmax < 1e6);
double wallstop = realdtmax + wallclock_time();
while( paras < nmax && (ww->nh == 0 || ww->wh[ww->nh-1]->get_t_max() < tmax) ) {
if(timelimit) {
double dt = wallstop - wallclock_time();
if(!hasdata(ins, dt))
break;
if(dt <= 0)
break;
}
if(ww->nh >= ww->maxnh) {
ww->maxnh = (ww->nh > ww->maxnh*2) ? ww->nh + 1 : ww->maxnh*2;
ww->wh = RenewN( ww->wh, worldbundleptr, ww->maxnh );
}
worldbundle *wb;
wb = read_bundle(*ins, ww->readflags & KIRA_VERBOSE);
if(wb == NULL) {
ins->close();
if(paras == 0) paras = -1;
break;
}
paras++;
}
if(paras > 0) {
ww->tmin = ww->wh[0]->get_t_min();
ww->tmax = ww->wh[ww->nh-1]->get_t_max();
}
return paras;
}
int kira_get_trange( struct dyndata *dd, struct stuff *st, double *tminp, double *tmaxp )
struct worldstuff *ww = (struct worldstuff *)dd->data;
if(ww == NULL) {
*tminp = *tmaxp = 0;
return 0;
}
*tminp = ww->tmin;
*tmaxp = ww->tmax;
return 1;
}
int kira_open( struct dyndata *dd, struct stuff *st, char *filename, int readflags)
{
ifstream *ins = new ifstream(filename);
if (!ins || !*ins) {
msg("Kira/starlab data file %s not found.", filename);
delete ins;
struct worldstuff *ww = NewN( worldstuff, 1 ); // from shmem -- avoid "new"
ww->maxnh = 1024;
ww->wh = NewN( worldbundleptr, ww->maxnh );
ww->ins = ins;
ww->tmin = ww->tmax = 0;
ww->tcur = ww->treq = -1.0; /* invalidate */
memset( dd, 0, sizeof(*dd) );
dd->data = ww;
dd->getspecks = kira_get_parti;
dd->trange = kira_get_trange;
dd->ctlcmd = kira_parse_args;
dd->draw = kira_draw;
dd->help = kira_help;
dd->free = NULL; /* XXX create a 'free' function someday */
dd->enabled = 1;
ww->nh = 0;
if(! (readflags & KIRA_READLATER)) {
kira_read_more( dd, st, ww->maxnh, 1e38, 1e38 );
if (ww->nh == 0) {
msg("Data file %s not in worldbundle format.", filename);
return 0;
}
}
ww->ih = 0;
ww->sl = NULL;
// If data field names aren't already assigned, set them now.
for(int i = 0; fieldnames[i] != NULL; i++) {
if(st->vdesc[st->curdata][i].name[0] == '\0')
strcpy(st->vdesc[st->curdata][i].name, fieldnames[i]);
}
// Register for H-R diagram display
static float xrange[] = { 5, 3 }; // initial range of log(T)
static float yrange[] = { -4, 6 }; // initial range of log(L)
static char xtitle[] = "log T";
static char ytitle[] = "log L";
plot = parti_register_plot( st, (void (*)(Fl_Plot*,void*,void*))kira_HRplot, dd );
if(plot) {
if(plot->x0() == plot->x1())
plot->xrange( xrange[0], xrange[1] );
if(plot->y0() == plot->y1())
plot->yrange( yrange[0], yrange[1] );
if(plot->xtitle() == NULL) plot->xtitle( xtitle );
if(plot->ytitle() == NULL) plot->ytitle( ytitle );
double get_parti_time( struct dyndata *dd, struct stuff *st ) {
return dd->data ? ((struct worldstuff *)(dd->data))->tcur : 0.0;
void set_parti_time( struct dyndata *dd, struct stuff *st, double reqtime ) {
if(dd->data == NULL)
static float log10_of( float v ) {
return v <= 0 ? 0 : log10f( v );
}
static float netTL( pdyn *b, float *totlum ) { // temp-weighted lum
if(b == NULL) return 0;
if(b->is_leaf()) {
*totlum += b->get_luminosity();
return b->get_temperature() * b->get_luminosity();
}
float totTL = 0;
for_all_daughters(pdyn, b, bb)
totTL += netTL( bb, totlum );
return totTL;
}
speck *add_speck(pdyn *b, pdyn *top, int addr, speck *sp, specklist *sl, worldstuff *ww)
sp->p.x[0] = b->get_pos()[0];
sp->p.x[1] = b->get_pos()[1];
sp->p.x[2] = b->get_pos()[2];
int id = b->get_index();
if(b->is_leaf()) {
if(id < ww->nleafsel) {
ww->unionsel |= ww->leafsel[id];
ww->interactsel = ww->intdest.wanton;
}
sp->val[SPECK_ID] = id;
sp->val[SPECK_TLOG] = log10_of( b->get_temperature() );
sp->val[SPECK_LUM] = b->get_luminosity();
id = -b->get_worldline_index();
if(id >= 0 || id + ww->maxstars < ww->maxleaves) {
static int yikes = 1;
if(yikes++ % 256 == 0) printf("OOPS: worldline_index %d not in range 1..%d\n",
-id, 2*ww->maxleaves);
int slotno = ww->maxstars + id;
ww->unionsel |= ww->leafsel[slotno];
sp->val[SPECK_ID] = id; // distinguish CM nodes
float totL = 0, totTL;
totTL = netTL( b, &totL );
sp->val[SPECK_LUM] = totL;
sp->val[SPECK_TLOG] = log10_of( totL == 0 ? 0 : totTL / totL );
sp->val[SPECK_MASSLOG] = log10_of( b->get_mass() * ww->massscale );
sp->val[SPECK_STYPE] = b->get_stellar_type(); /* see starlab/inc/star_support.h */
sp->val[SPECK_NCLUMP] = 0; // complete (add n_leaves) in kira_to_parti
// DON'T offset unperturbed binaries by 100 in "nclump" field.
//if (b != top) {
// pdyn *bb = b;
// while (bb->get_elder_sister())
// bb = bb->get_elder_sister();
// if (bb->get_kepler()) sp->val[SPECK_NCLUMP] += 100;
//}
sp->val[SPECK_ROOTID] = (top->is_leaf())
? top->get_index()
: -top->get_worldline_index();
// Addr is 0 for top-level leaf, 1 for top-level CM,
// constructed from parent addr otherwise.
if(id != 0 && ww->tracking == id) {
Point newpos;
newpos.x[0] = b->get_pos()[0];
newpos.x[1] = b->get_pos()[1];
newpos.x[2] = b->get_pos()[2];
if(ww->wastracking) {
Point delta;
vsub( &delta, &newpos,&ww->trackpos );
parti_nudge_camera( &delta );
} else {
kira_track_break( ww, &newpos );
}
ww->trackpos = newpos;
ww->wastracking = 1;
}
speck *assign_specks(pdyn *b, pdyn *top, int addr, speck *sp, specklist *sl, worldstuff *ww)
vector oldpos = b->get_pos();
if(b->is_leaf()
|| ww->treenodes == KIRA_ON
|| (ww->treenodes == KIRA_ROOTS && b == top)) {
sp = add_speck(b, top, addr, sp, sl, ww);
if(b->is_leaf())
ww->leafcount++;
}
if(!b->is_leaf() &&
|| (ww->treerings == KIRA_ROOTS && b == top))) {
struct speck *marksp = ww->marksp;
ww->marksp = add_speck(b, top, addr, marksp, ww->marksl, ww);
// find our two children
pdyn *b1 = b->get_oldest_daughter();
pdyn *b2 = b1 ? b1->get_younger_sister() : NULL;
if(b2) {
vector sep = b1->get_pos() - b2->get_pos();
real dist = sqrt(sep * sep);
break;
case KIRA_RINGA:
{
vector vel = b1->get_vel() - b2->get_vel();
real speed2 = vel * vel;
real E = 0.5 * speed2 - mass / dist;
size = -mass / (2*E); // semimajor axis
}
}
marksp->val[SPECK_RINGSIZE] = size;
marksp->val[SPECK_MU] = b1->get_mass() / mass;
marksp->val[SPECK_SEPVEC] = sep[0];
marksp->val[SPECK_SEPVEC+1] = sep[1];
marksp->val[SPECK_SEPVEC+2] = sep[2];
if(b != top)
b->set_pos( oldpos );
struct specklist *kira_to_parti(pdyn *root, struct dyndata *dd, struct stuff *st, struct specklist *sl)
worldstuff *ww = (worldstuff *)dd->data;
// Use NewN to allocate these in (potentially) shared memory,
// rather than new which will allocate in local malloc pool.
// Needed for virdir/cave version.
int w, wmax = 0;
for(int k = 0; k < ww->nh; k++) {
if(ww->wh[k] && wmax < (w = ww->wh[k]->get_nw()))
wmax = w;
}
ww->maxstars = smax + wmax + smax; // leave a bit of room for growth
ww->maxleaves = smax;
sl->bytesperspeck = SMALLSPECKSIZE( SPECK_NDATAFIELDS );
sl->specks = NewNSpeck(sl, ww->maxstars);
sl->sel = NewN( SelMask, ww->maxstars );
sl->nsel = ww->maxstars;
memset(sl->sel, 0, ww->maxstars*sizeof(int));
marksl->sel = NewN( SelMask, ww->maxmarks );
marksl->nsel = ww->maxmarks;
marksl->specks = NewNSpeck(marksl, ww->maxmarks);
marksl->special = MARKERS;
sl->next = marksl;
marksl->next = NULL;
ww->leafsel = NewN( SelMask, ww->maxstars );
ww->nleafsel = ww->maxstars;
memset(ww->leafsel, 0, ww->maxstars*sizeof(SelMask));
ww->trails = NewN( trailhead, ww->maxstars );
memset(ww->trails, 0, ww->maxstars*sizeof(trailhead));
selinit( &ww->intdest );
selinit( &ww->intsrc );
specks_all_picks( st, kira_picked, dd );
if(ww->myselseq < sl->selseq && sl->sel != NULL) {
// somebody changed something in the global select bits, so
// gather all stars' sel-bits back into our local copy.
struct speck *tsp = sl->specks;
for(int k = 0; k < sl->nsel; k++) {
int id = (int)tsp->val[SPECK_ID];
if(id > 0 && id < ww->nleafsel) {
ww->leafsel[id] = sl->sel[k];
} else if(id < 0 && id + ww->maxstars >= ww->maxleaves) {
ww->leafsel[ww->maxstars + id] = sl->sel[k];
}
tsp = NextSpeck(tsp, sl, 1);
}
sl->colorseq = sl->sizeseq = sl->threshseq = 0;
marksl->colorseq = marksl->sizeseq = marksl->threshseq = 0;
int i;
struct vald *vd = ww->vd;
for(i = 0, vd = ww->vd; i < SPECK_NDATAFIELDS; i++, vd++) {
vd->min = 1e9;
vd->max = -1e9;
vd->sum = 0;
}
int ntotal = 0;
ww->interactsel = ww->unionsel = 0;
int k;
int nleaves = ww->leafcount - nl;
int count = sl->nspecks - ns;
ww->unionsel |= ww->interactsel;
/* For all nodes in this subtree, ... */
int id = (int)tsp->val[SPECK_ID];
if(id < 0) {
/* For CM nodes, negate nclump value. */
tsp->val[SPECK_NCLUMP] = -tsp->val[SPECK_NCLUMP];
/* Also, propagate all leaves' set membership to CM nodes */
sl->sel[ntotal] = ww->unionsel;
} else if(ww->interactsel && id < ww->nleafsel) {
/* For leaf nodes, if at least one star in this
* group is in the interaction set,
* then add all other group members to it.
*/
sl->sel[ntotal] = ww->leafsel[id] |= ww->interactsel;
/* Making trails? */
}
if(ww->trailsel.use != SEL_NONE &&
SELECTED( ww->leafsel[id], &ww->trailsel )) {
kira_add_trail( st, ww, id, tsp );
} else if(ww->trailonly) {
kira_erase_trail( st, ww, id );
for(i = 0, vd = ww->vd; i <= SPECK_STYPE; i++, vd++) {
float v = tsp->val[i];
if(vd->min > v) vd->min = v;
if(vd->max < v) vd->max = v;
vd->sum += v;
}
ntotal++;
/* For all marks (rings, etc.) in this subtree */
count = marksl->nspecks - mns;
for(k = 0; k < count; k++) {
msp->val[SPECK_NCLUMP] += nleaves;
if(msp->val[0] < 0)
msp->val[SPECK_NCLUMP] = -msp->val[SPECK_NCLUMP]; // negate nclump for CM nodes
marksl->sel[mns+k] = ww->unionsel;
if(ww->wastracking < 0) ww->wastracking = 0; // Detach if tracked pcle not found now
for(i = 0, vd = ww->vd; i < SPECK_NDATAFIELDS && i < MAXVAL; i++, vd++) {
struct valdesc *tvd = &st->vdesc[ st->curdata ][ i ];
tvd->nsamples = ntotal;
if(vd->min > vd->max) {
vd->min = vd->max = 0;
tvd->nsamples = 0;
}
tvd->min = vd->min;
tvd->max = vd->max;
tvd->sum = vd->sum;
tvd->mean = tvd->nsamples > 0 ? vd->sum / tvd->nsamples : 0;
}
specks_reupdate( st, st->sl );
struct specklist *kira_get_parti( struct dyndata *dd, struct stuff *st, double realtime )
struct worldstuff *ww = (struct worldstuff *)dd->data;
return NULL;
ww->treq = realtime;
if (realtime < ww->tmin) realtime = ww->tmin;
if (realtime > ww->tmax) realtime = ww->tmax;
return ww->sl;
int ih = ww->ih;
int nh = ww->nh;
worldbundle *wb = ww->wh[ih];
for(; realtime > wb->get_t_max() && ih < nh-1; ih++, wb = ww->wh[ih])
;
for(; realtime < wb->get_t_min() && ih > 0; ih--, wb = ww->wh[ih])
;
#if OLDTREE
pdyn *root = create_interpolated_tree(wb, realtime);
ww->sl = kira_to_parti(root, dd, st, ww->sl);
rmtree(root);
#else
pdyn *root = create_interpolated_tree2(wb, realtime);
ww->sl = kira_to_parti(root, dd, st, ww->sl);
#endif
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
void kira_add_trail( struct stuff *st, worldstuff *ww, int id, struct speck *sp )
{
if(id < 0) id = ww->maxstars + id;
if(id <= 0 || id > ww->maxstars || ww->maxtrail <= 0) return;
struct trailhead *th = &ww->trails[id];
int bps = ww->sl->bytesperspeck;
if(th->maxtrail <= 0) {
th->maxtrail = ww->maxtrail;
th->specks = (struct speck *)NewN( char, th->maxtrail*bps );
th->ntrails = 0;
th->next = 0;
}
if(th->ntrails > th->maxtrail)
th->ntrails = th->maxtrail;
if(th->next < 0 || th->next >= th->maxtrail)
th->next = 0;
struct speck *tsp = (struct speck *)(((char *)th->specks) + th->next*bps);
memcpy( tsp, sp, bps );
vsub( &tsp->p, &tsp->p, &ww->trackpos );
if(++th->next >= th->maxtrail)
th->next = 0;
if(++th->ntrails > th->maxtrail)
th->ntrails = th->maxtrail;
if(id >= ww->maxtrailno)
ww->maxtrailno = id+1;
}
void kira_erase_trail( struct stuff *st, worldstuff *ww, int id )
{
if(id < 0) id = ww->maxstars + id;
if(id <= 0 || id > ww->maxstars || ww->maxtrail <= 0) return;
struct trailhead *th = &ww->trails[id];
th->next = th->ntrails = 0;
if(id+1 == ww->maxtrailno) {
while(ww->maxtrailno > 0 && ww->trails[ww->maxtrailno-1].ntrails == 0)
ww->maxtrailno--;
}
}
void kira_track_break( struct worldstuff *ww, Point *newtrack )
{
Point incr;
vsub( &incr, newtrack, &ww->trackpos );
for(int id = 0; id < ww->maxstars; id++) {
struct trailhead *th = &ww->trails[id];
if(th->ntrails > 0) {
struct speck *sp = th->specks;
for(int i = 0; i < th->maxtrail; i++) {
vsub( &sp->p, &sp->p, &incr );
sp = NextSpeck(sp, ww->sl, 1);
}
}
}
// ww->trackpos = *newtrack; no, let caller do that.
}
void kira_maxtrail( struct dyndata *dd, struct stuff *st, int newmax )
{
int i;
struct worldstuff *ww = (struct worldstuff *)dd->data;
ww->maxtrail = newmax;
if(ww->trails == NULL) return;
int bps = ww->sl->bytesperspeck;
char *spare = (char *)malloc( bps * newmax );
for(i = 0; i < ww->maxstars; i++) {
struct trailhead *th = &ww->trails[i];
if(th->maxtrail == 0)
continue;
if(th->ntrails > newmax) th->ntrails = newmax;
int first = (th->next + th->maxtrail - th->ntrails) % th->maxtrail;
char *base = (char *)th->specks;
if(first < th->next) {
memmove( spare, base + first*bps, th->ntrails*bps );
} else {
/* rearrange two pieces: 0..next-1 first..max-1
* into 0..keep-1 keep..keep+next-1
*/
int keep = th->maxtrail - first;
memcpy( spare, base + first*bps, keep*bps );
memmove( spare + keep*bps, base, th->next*bps );
}
Free( th->specks );
th->specks = (struct speck *)NewN( char, newmax*bps );
memcpy( th->specks, spare, th->ntrails*bps );
th->next = th->ntrails;
th->maxtrail = newmax;
}
free(spare);
}
int kira_draw( struct dyndata *dd, struct stuff *st, struct specklist *slhead, Matrix *Tc2w, float radperpix )
struct worldstuff *ww = (struct worldstuff *)dd->data;
if(ww == NULL || !dd->enabled)
return 0;
int inpick = st->inpick;
int slno;
static Point zero = {0,0,0};
Point fwdvec = {0,0,-radperpix}; // scaled by pixels-per-radian
float fwdd;
vtfmpoint( &eyepoint, &zero, Tc2w );
vtfmvector( &fwd, &fwdvec, Tc2w );
fwdd = -vdot( &eyepoint, &fwd );
#define MAXRING 32
Point fan[MAXRING];
int nring = 24;
int i;
for(i = 0; i < nring; i++) {
float th = 2*M_PI*i / nring;
vcomb( &fan[i], cos(th),(Point *)&Tc2w->m[0*4+0],
sin(th),(Point *)&Tc2w->m[1*4+0] );
}
int specks_slno = 0;
for(sl = slhead, slno = 1; sl != NULL; sl = sl->next, slno++) {
if(sl == ww->sl) {
specks_slno = slno;
} else if(sl->special == MARKERS) {
if(inpick) {
glLoadName(slno);
glPushName(0);
}
for(i = 0; i < ns; i++, sp = NextSpeck(sp, sl, 1)) {
float unitperpix = vdot( &fwd, &sp->p ) + fwdd;
if(unitperpix <= 0) continue;
if( !SELECTED(sl->sel[i], &st->seesel) )
continue;
if(ww->treerings == KIRA_ON ||
(ww->treerings == KIRA_ROOTS && sp->val[SPECK_TREEADDR] == 1)) {
float ringpixels = fabs(sp->val[SPECK_RINGSIZE] * ww->ringscale) / unitperpix;
if(ringpixels < ww->ringmin) ringpixels = ww->ringmin;
if(ringpixels > ww->ringmax) ringpixels = ww->ringmax;
float rring = ringpixels * unitperpix;
int step = ringpixels>20 ? 1 : ringpixels>8 ? 2 : 3;
if(inpick) glLoadName(i);
else glColor3ubv( (GLubyte *)&sp->rgba );
glBegin( GL_LINE_LOOP );
for(int k = 0; k < nring; k+=step)
glVertex3f( sp->p.x[0] + rring*fan[k].x[0],
sp->p.x[1] + rring*fan[k].x[1],
sp->p.x[2] + rring*fan[k].x[2] );
glEnd();
}
if(ww->treearcs != KIRA_OFF) {
float mu = sp->val[SPECK_MU];
float *sep = &sp->val[SPECK_SEPVEC];
if(inpick) glLoadName(i);
else glColor3ubv( (GLubyte *)&sp->rgba );
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
glBegin( GL_LINES );
if(ww->treearcs != KIRA_TICK) {
glVertex3f(
sp->p.x[0] - mu*sep[0],
sp->p.x[1] - mu*sep[1],
sp->p.x[2] - mu*sep[2] );
glVertex3f(
sp->p.x[0] + (1-mu)*sep[0],
sp->p.x[1] + (1-mu)*sep[1],
sp->p.x[2] + (1-mu)*sep[2] );
}
if(ww->treearcs == KIRA_CROSS || ww->treearcs == KIRA_TICK) {
Point cr;
vcross( &cr, &unitfwd, (Point *)sep );
float sepsep = vdot((Point *)sep, (Point *)sep);
float sepfwd = vdot((Point *)sep, &unitfwd);
float scaleby = halftickscale / sqrt(1 - sepfwd*sepfwd/sepsep);
glVertex3f(
sp->p.x[0] - scaleby*cr.x[0],
sp->p.x[1] - scaleby*cr.x[1],
sp->p.x[2] - scaleby*cr.x[2] );
glVertex3f(
sp->p.x[0] + scaleby*cr.x[0],
sp->p.x[1] + scaleby*cr.x[1],
sp->p.x[2] + scaleby*cr.x[2] );
}
glEnd();
}
if(inpick)
glPopName();
glLineWidth( ww->trailpsize );
glEnable( GL_BLEND );
glBlendFunc( GL_SRC_ALPHA, GL_ONE );
int alpha = 0;
((GLubyte *)&alpha)[3] = (int) (255 * ww->trailalpha);
if(inpick) {
glPushName( specks_slno );
glPushName( 0 );
}
for(i = 0; i < ww->maxtrailno; i++) {
struct trailhead *th = &ww->trails[i];
if(th->ntrails == 0 || th->maxtrail <= 0)
continue;
if(inpick)