From: Gordon L. K. <gl...@uc...> - 2009-12-01 16:52:47
|
Hello, I wanted to pick up this thread from a little while ago, in particular the stuff about gageProbe'ing many values at once. I have a planned API for how to do this, and wanted feedback from others. On Apr 23, 2008, at 2:21 PM, Carlos Scheidegger wrote: >> hi- >> >>> I needed to extend gage to probe at many points in a single >>> function >>> call. I'm writing a Python wrapper for teem (so far, nrrd and gage >>> are >>> essentially wrapped). >> >> sweet. Jorik Blaas (hi Jorik!) has also done some Teem wrapping >> work ... > > For whatever's worth, I'm using ctypes and a mostly manual > translation of the headers ... > >>> The issue is that Python has a huge function- >>> call overhead, so calling gageProbe from python many times is a bad >>> idea if I can avoid it. This led me to write gageProbeMany, a simple >>> function to perform a set of gageProbe queries: >>> >>> GAGE_EXPORT int gageProbeMany(gageContext *ctx, >>> Nrrd *nout[GAGE_PERVOLUME_MAXNUM][GAGE_ITEM_MAX], >>> const Nrrd *nin); The concern I voiced at the time is that this entrenched GAGE_PERVOLUME_MAXNUM and GAGE_ITEM_MAX in the API, which I didn't want to do, because in the long run neither of these constants should exist. Indeed, GAGE_PERVOLUME_MAXNUM has since been removed (and will stay removed for the 1.11 release), because the new scale-space functionality in Gage allows more per-volumes than had ever been anticipated before. Code that relies on GAGE_PERVOLUME_MAXNUM won't compile any more, nor should it, because it probably won't work anyway. I'll buy lunch for anyone for whom this is actually a problem (hey Carlos); the documented usgae of gage never refers to this. GAGE_ITEM_MAX may also eventually go away. Below is my proposed addition to gage.h, with a little usage example. Feedback welcome. Gordon /* ******** gageMultiItem ** ** used (along with gageMultiQuery) to represent a list of items that ** we want to learn from a single pervolume, at every point given in a ** Nrrd of probe positions. The answers for those items are ** concatenated along axis 0 of the single output Nrrd "nans", which ** for simplicity of access is stored within here as well. */ typedef struct { unsigned int itemNum; /* number of items we want to learn, allocated length of item[], ansDirect[], and ansLen[] */ int *item; /* all items for this pervolume */ double *ansDir; /* ansDir[i] points to answer i */ unsigned int *ansLen; /* ansLen[i] length of answer i */ /* ========== OUTPUT */ Nrrd *nans; /* output Nrrd */ } gageMultiItem; /* ******** gageMultiQuery ** ** used to represent queries, and the organization of those queries ** that should be probed at a given Nrrd of probe positions. The ** basic idea is that we need to be able to represent *1* the goal of ** having two or more answers concatenated together (along axis 0), ** in a single Nrrd (of the sort that one might pass to nrrdApply2DLut) ** AND *2* the goal of having different answers stored in different Nrrds. ** This is handled by having, for each pervolume, an array (possibly ** length 1) of queries, and generating one Nrrd per query. */ typdedef struct { unsigned int pvlNum, /* number of perVolumes in the context for which we represent the queries; this is also the allocated length of queryNum[ii] and query[ii] */ *queryNum; /* queryNum[ii] is the number of queries for pvl[ii]; which is also (due to use of gageMultiItem) the number of answer Nrrds for pvl[ii] */ gageMultiItem ***query; /* query[ii] is array of gageMultiItem*s for pvl[ii], each is query[ii][nn]; nn = 0..queryNum[ii] */ } gageMultiQuery; gageMultiItem *gageMultiItemNew(); int gageMultiItemSet_va(gageMultiItem *gmi, unsigned int itemNum, ... /* itemNum items */); gageMultiItem *gageMultiItemNix(gageMultiItem *gmi); gageMultiItem *gageMultiItemNuke(gageMultiItem *gmi); gageMultiQuery *gageMultiQueryNew(gageContext *gctx); int gageMultiQueryAdd_va(gageMultiQuery *gmq, unsigned int pvlIdx, unsigned int queryNum, ... /* queryNum gageMultiItem* */); int gageMultiProbe(gageMultiQuery *gmq, gageContext *gctx, const Nrrd *npos); int gageMultiProbeSpace(gageMultiQuery *gmq, gageContext *gctx, const Nrrd *npos, int indexSpace, int clamp); gageMultiQuery *gageMultiQueryNix(gageMultiQuery *gmq); gageMultiQuery *gageMultiQueryNuke(gageMultiQuery *gmq); example usage: ... set a tensor pervolume in gageContext gctx ... gageMultiItem *gmi0, *gmi1; gmi0 = gageMultiItemNew(); gmi1 = gageMultiItemNew(); gageMultiItemSet_va(gmi0, 2, tenGageFA, tenGageMode); gageMultiItemSet_va(gmi1, 1, tenGageEvec0); gageMultiQuery *gmq; gmq = gageMuliQueryNew(gctx); gageMultiQueryAdd_va(gmq, 0, 2, gmi0, gm1); /* a second pervolume would be handled with gageMultiQueryAdd_va(gmq, 1, n, ...) */ gageMultiProbe(gmq, gctx, npos); nrrdSave("fa-mode.nrrd", gmi0->nans, NULL); nrrdSave("evec0.nrrd", gmi1->nans, NULL); gageMultiQueryNuke(gmq); |