StEventUtilities/StEventHelper.cxx

00001 #include <stdio.h> 00002 #include <stdlib.h> 00003 #include <assert.h> 00004 #include "TROOT.h" 00005 #include "TClass.h" 00006 #if ROOT_VERSION_CODE < 331013 00007 #include "TCL.h" 00008 #else 00009 #include "TCernLib.h" 00010 #endif 00011 #include "TMath.h" 00012 #include "TBaseClass.h" 00013 #include "TDataMember.h" 00014 #include "TMethod.h" 00015 #include "TMethodArg.h" 00016 #include "TDataType.h" 00017 #include "Api.h" 00018 #include "TMemberInspector.h" 00019 #include "TExMap.h" 00020 #include "TCollection.h" 00021 #include "TRegexp.h" 00022 #include "TRandom.h" 00023 #include "TError.h" 00024 #include "TPoints3DABC.h" 00025 #include "TSystem.h" 00026 #include "TPad.h" 00027 #include "TView.h" 00028 00029 #include "StEvent.h" 00030 #include "StObject.h" 00031 #include "StHit.h" 00032 #include "StTrack.h" 00033 #include "StTrackNode.h" 00034 #include "StVertex.h" 00035 #include "StTrackGeometry.h" 00036 #include "StTrackDetectorInfo.h" 00037 // For L3 filter 00038 #include "StarClassLibrary/BetheBloch.h" 00039 #include "StEvent/StDedxPidTraits.h" 00040 // For coloring filter 00041 #include "StEventTypes.h" 00042 #include "StProbPidTraits.h" 00043 #include "StTpcDedxPidAlgorithm.h" 00044 #include "THelixTrack.h" 00045 00046 #define __EVENTHELPER_ONLY__ 00047 #include "StEventHelper.h" 00048 #undef __EVENTHELPER_ONLY__ 00049 #include <map> 00050 00051 #include "StSvtHitCollection.h" 00052 #include "StSvtBarrelHitCollection.h" 00053 #include "StSvtLadderHitCollection.h" 00054 #include "StSvtWaferHitCollection.h" 00055 00056 void Break(){printf("InBreak\n");} 00057 00058 std::map<long,long> myMap; 00059 typedef std::pair <long,long> MyPair; 00060 std::map <long,long> :: const_iterator myFinder; 00061 00062 00063 00064 00065 class StEventInspector : public TMemberInspector { 00066 public: 00067 StEventInspector(TExMap *map,Int_t &count,const char *opt=""); 00068 virtual ~StEventInspector(){delete fSkip;}; 00069 virtual void Inspect(TClass* cl, const char* parent, const char* name, const void* addr); 00070 void CheckIn(TObject *obj,const char *bwname=""); 00071 00072 Int_t &fCount; 00073 TExMap *fMap; 00074 TRegexp *fSkip; 00075 TString fOpt; 00076 00077 }; 00078 00079 //______________________________________________________________________________ 00080 StEventInspector::StEventInspector(TExMap *map,Int_t &count,const char *opt):fCount(count) 00081 { 00082 fMap = map; 00083 fSkip = 0; 00084 fOpt = opt; 00085 if (fOpt.Length()) fSkip = new TRegexp(fOpt.Data()); 00086 } 00087 //______________________________________________________________________________ 00088 void StEventInspector::Inspect(TClass* kl, const char* tit , const char* name, const void* addr) 00089 { 00090 if(tit && strchr(tit,'.')) return ; 00091 00092 TString ts; 00093 00094 if (!kl) return; 00095 if (name[0] == '*') name++; 00096 int ln = strcspn(name,"[ "); 00097 TString iname(name,ln); 00098 const char *iName=iname.Data(); 00099 if (iName[1]=='P' && strcmp(iName,"fParent" )==0) return; 00100 if (iName[0]=='G' && strcmp(iName,"G__virtualinfo")==0) return; 00101 00102 G__ClassInfo *classInfo = (G__ClassInfo *)kl->GetClassInfo(); 00103 if (!classInfo) return; 00104 G__ClassInfo &cl = *classInfo; 00105 00106 00107 // Browse data members 00108 G__DataMemberInfo m(cl); 00109 int found=0; 00110 const char *mName=0; 00111 while (m.Next()) { // MemberLoop 00112 mName = m.Name(); 00113 if (mName[1] != iName[1]) continue; 00114 if (strcmp(mName,iName) ) continue; 00115 found = 1; break; 00116 } 00117 assert(found); 00118 00119 // we skip: non TObjects 00120 // - the member G__virtualinfo inserted by the CINT RTTI system 00121 00122 long prop = m.Property() | m.Type()->Property(); 00123 if (prop & G__BIT_ISFUNDAMENTAL) return; 00124 if (prop & G__BIT_ISSTATIC) return; 00125 if (prop & G__BIT_ISENUM) return; 00126 if (strcmp(m.Type()->Fullname(),"TObject") && !m.Type()->IsBase("TObject")) 00127 return; 00128 00129 int size = sizeof(void*); 00130 if (!(prop&G__BIT_ISPOINTER)) size = m.Type()->Size(); 00131 00132 int nmax = 1; 00133 if (prop & G__BIT_ISARRAY) { 00134 for (int dim = 0; dim < m.ArrayDim(); dim++) nmax *= m.MaxIndex(dim); 00135 } 00136 00137 for(int i=0; i<nmax; i++) { 00138 char *ptr = (char*)addr + i*size; 00139 TObject *obj = (prop&G__BIT_ISPOINTER) ? *((TObject**)ptr) : (TObject*)ptr; 00140 if (!obj) continue; 00141 const char *bwname = obj->ClassName(); 00142 if (!bwname[0] || strcmp(bwname,obj->ClassName())==0) { 00143 bwname = name; 00144 int l = strcspn(bwname,"[ "); 00145 if (bwname[l]=='[') { 00146 char cbuf[12]; sprintf(cbuf,"[%02d]",i); 00147 ts.Replace(0,999,bwname,l); 00148 ts += cbuf; 00149 bwname = (const char*)ts; 00150 } 00151 } 00152 00153 CheckIn(obj,bwname); 00154 00155 } 00156 00157 } 00158 //______________________________________________________________________________ 00159 void StEventInspector::CheckIn(TObject *obj,const char *bwname) 00160 { 00161 if (!obj) return; 00162 TObject *inobj=0; 00163 if (obj->InheritsFrom(StRefArray::Class())) return; 00164 if (obj->InheritsFrom( StObjLink::Class())) return; 00165 int n; 00166 if (fSkip && (fSkip->Index(obj->ClassName(),&n)>=0)) return; 00167 00168 if (obj->InheritsFrom(TCollection::Class())){ 00169 TCollection *tcol = (TCollection*)obj; 00170 TIter next(tcol); 00171 while ((inobj=next())) {CheckIn(inobj);} 00172 return; 00173 } 00174 00175 if (obj->InheritsFrom(StXRef::Class())){ 00176 00177 Long_t &inmap = (*fMap)(TMath::Hash(&obj,sizeof(void*)),(Long_t)obj); 00178 00179 myFinder = myMap.find((long)obj); 00180 assert((inmap==0) == (myFinder==myMap.end())); 00181 00182 if (inmap) return; 00183 myMap.insert(MyPair((long)obj,1)); 00184 inmap = 1;fCount++; 00185 } 00186 00187 if (obj->InheritsFrom(StStrArray::Class())){ 00188 // if (obj->IsA()==StSPtrVecTrackNode::Class()) Break(); 00189 if (((StStrArray*)obj)->size()) { 00190 00191 Long_t &inmap = (*fMap)(TMath::Hash(&obj,sizeof(void*)),(Long_t)obj); 00192 00193 myFinder = myMap.find((long)obj); 00194 assert((inmap==0) == (myFinder==myMap.end())); 00195 if (inmap) return; 00196 myMap.insert(MyPair((long)obj,2)); 00197 00198 inmap = 2; fCount++; 00199 int vecobj = ( obj->IsA() == StSPtrVecObject::Class()); 00200 // printf("%8d %p %s::%s\n",fLevel,(void*)obj,obj->GetName(),bwname); 00201 StStrArray *arr = (StStrArray*)obj; 00202 int sz = arr->size(); 00203 for (int idx=0;idx<sz; idx++) { 00204 inobj = arr->at(idx); 00205 Int_t count = fCount; 00206 CheckIn(inobj); 00207 if (count==fCount && !vecobj) break; //if no action was made, no sense to continue 00208 } 00209 return; 00210 } } 00211 00212 char cbuf[1000];*cbuf=0; 00213 StEventInspector insp(fMap,fCount); 00214 obj->ShowMembers(insp,cbuf); 00215 } 00216 //______________________________________________________________________________ 00217 ClassImp(StEventHelper) 00218 //______________________________________________________________________________ 00219 StEventHelper::StEventHelper(const TObject *evt,const char *opt) 00220 { 00221 fMap = new TExMap(10000); 00222 myMap.clear(); 00223 fObject = 0; 00224 Reset(evt,opt); 00225 } 00226 //______________________________________________________________________________ 00227 StEventHelper::~StEventHelper() 00228 { 00229 myMap.clear(); 00230 delete fMap; fMap=0; 00231 Clear(); 00232 } 00233 //______________________________________________________________________________ 00234 void StEventHelper::Clear(Option_t *opt) 00235 { 00236 } 00237 //______________________________________________________________________________ 00238 void StEventHelper::Reset(const TObject *evt,const char *opt) 00239 { 00240 fObject = (TObject *)evt; 00241 Clear(); 00242 myMap.clear(); 00243 fMap->Delete(); 00244 if (!fObject) return; 00245 int kount=0; 00246 StEventInspector insp(fMap,kount,opt); 00247 char cbuf[1024]; 00248 00249 fObject->ShowMembers(insp,cbuf); 00250 } 00251 //______________________________________________________________________________ 00252 int StEventHelper::Kind(const TObject *to) 00253 { 00254 static TClass *klass=0; 00255 static int who=0; 00256 int kind = 0; 00257 TClass *myClass=to->IsA(); 00258 if (myClass!=klass) { 00259 klass = myClass; who = 0; 00260 if (klass->InheritsFrom( StHit::Class())) { who=kHIT;} 00261 else if (klass->InheritsFrom( StTrack::Class())) { who=kTRK;} 00262 else if (klass->InheritsFrom(StPtrVecHit::Class())) { who=kHRR;} 00263 else if (klass->InheritsFrom( StVertex::Class())) { who=kVTX;} 00264 else if (klass->InheritsFrom( TObjArray::Class())) { who=kTRR;} 00265 } 00266 kind = who; 00267 if (kind==kHIT) { 00268 StHitHelper hh((StHit*)to); 00269 if (hh.IsUsed()) {kind|=kUSE;} else {kind|=kUNU;} 00270 if (hh.IsFit ()) kind|=kFIT; 00271 return kind; 00272 } 00273 return kind; 00274 } 00275 //______________________________________________________________________________ 00276 THelixTrack *StEventHelper::MyHelix(THelixTrack *myHlx,const StHelixD *evHlx) 00277 { 00278 if (!myHlx) myHlx= new THelixTrack; 00279 double curv = evHlx->curvature(); 00280 double phase = evHlx->phase(); 00281 double dip = evHlx->dipAngle(); 00282 int h = evHlx->h(); 00283 00284 double myDir[3]; 00285 myDir[0]= -sin(phase)*cos(dip); 00286 myDir[1]= cos(phase)*cos(dip); 00287 myDir[2]= sin(dip); 00288 if (h<0) {myDir[0]=-myDir[0]; myDir[1]=-myDir[1];} 00289 double myX[3]; 00290 myX[0]= evHlx->x(0.); 00291 myX[1]= evHlx->y(0.); 00292 myX[2]= evHlx->z(0.); 00293 00294 myHlx->Set(myX,myDir,curv*h); 00295 return myHlx; 00296 } 00297 //______________________________________________________________________________ 00298 void StEventHelper::ls(Option_t* option) const 00299 { 00300 typedef struct { int nb; int sz; const char *tenant; } QWE; 00301 QWE *qwe=0; 00302 00303 TExMap map; 00304 TExMapIter it(fMap); 00305 Long_t key,val; 00306 while( it.Next(key,val) ) { 00307 if (val != 2) continue; 00308 StStrArray *a = (StStrArray *)(key); 00309 Long_t &cnt = map((Long_t)a->IsA()); 00310 // printf("%s %p\n",a->ClassName(),(void*)a); 00311 if (!cnt) { 00312 qwe = new QWE; 00313 cnt = (Long_t)qwe; 00314 qwe->nb=0; qwe->sz=0;qwe->tenant=0; 00315 } 00316 qwe = (QWE*)cnt; 00317 qwe->nb++; qwe->sz += a->size(); 00318 if (qwe->tenant==0 && a->size()) { 00319 TObject *to = a->front(); 00320 if (to) qwe->tenant = to->ClassName(); 00321 } 00322 00323 } 00324 TExMapIter itt(&map); 00325 printf("\n StEvent(%p)\n",(void*)fObject); 00326 00327 while( itt.Next(key,val) ) { 00328 TObject *kl = (TObject *)key; 00329 qwe = (QWE*)val; 00330 printf ("%8d(%8d) - %s (%s)\n",qwe->nb,qwe->sz,kl->GetName(),qwe->tenant); 00331 delete qwe; 00332 } 00333 printf("\n"); 00334 00335 } 00336 //______________________________________________________________________________ 00337 TObjArray *StEventHelper::SelConts(const char *sel) 00338 { 00339 TObjArray *tarr = new TObjArray; 00340 TRegexp reg(sel); 00341 00342 TExMapIter it(fMap); 00343 Long_t key,val; 00344 while( it.Next(key,val) ) { 00345 if (val == 1) continue; 00346 StStrArray *a = (StStrArray *)(key); 00347 if(a->size()==0) continue; 00348 int n =0; 00349 if (reg.Index(a->ClassName(),&n)<0) continue; 00350 tarr->Add(a); 00351 } 00352 return tarr; 00353 } 00354 //______________________________________________________________________________ 00355 TObjArray *StEventHelper::SelTracks(const char*,int flag) 00356 { 00357 int trackTypes[]= {global, primary, tpt, secondary, estGlobal, estPrimary,-1}; 00358 00359 TObjArray *conts = SelConts("^StSPtrVecTrackNode$"); 00360 TObjArray *traks = new TObjArray(); 00361 Int_t ilast = conts->GetLast(); 00362 for (int idx=0;idx<=ilast;idx++) { 00363 StObjArray *arr = (StObjArray *)conts->At(idx); 00364 if (!arr) continue; 00365 int ntrk = arr->size(); 00366 if (!ntrk) continue; 00367 for (int itrk=0;itrk<ntrk;itrk++) { 00368 StTrackNode *tn = (StTrackNode*)arr->at(itrk); 00369 if (!tn) continue; 00370 StTrack *trk = 0;int ity; int bty=kTGB; 00371 for (int jty=0;(ity=trackTypes[jty])>=0;jty++,bty<<=1){ 00372 if (!(flag&bty)) continue; 00373 trk=tn->track(ity); 00374 if (!trk) continue; 00375 if (trk->IsZombie()) continue; 00376 // See: StRoot/St_base/StObject.h also 00377 if ((flag&kMark2Draw) && !trk->TestBit(kMark2Draw)) continue; 00378 traks->Add(trk); 00379 }//end track types 00380 }//end StTrackNode's 00381 }// end StSPtrVecTrackNode's 00382 delete conts; 00383 return traks; 00384 } 00385 //______________________________________________________________________________ 00386 TObjArray *StEventHelper::SelHits(const char *RegEx, Int_t un, Int_t flag) 00387 { 00388 // un == used +2*nonused 00389 00390 TObjArray *conts = SelConts(RegEx); 00391 TObjArray *hits = new TObjArray(); 00392 Int_t ilast = conts->GetLast(); 00393 00394 for (int idx=0;idx<=ilast;idx++) { 00395 StObjArray *arr = (StObjArray *)conts->At(idx); 00396 if (!arr) continue; 00397 int sz = arr->size(); 00398 if (!sz) continue; 00399 if (!arr->at(0)->InheritsFrom(StHit::Class())) continue; 00400 for(int ih=0;ih<sz; ih++) { 00401 StHit *hit = (StHit*)arr->at(ih); 00402 if (!hit) continue; 00403 if (hit->IsZombie()) continue; 00404 // See: StRoot/St_base/StObject.h also 00405 if ((flag&kMark2Draw) && !hit->TestBit(kMark2Draw)) continue; 00406 int used = (hit->trackReferenceCount()!=0); 00407 int take = 0; 00408 if ( used && (un&kUSE)) take++; 00409 if (!used && (un&kUNU)) take++; 00410 if (take) hits->Add(hit); 00411 } 00412 } 00413 delete conts; 00414 return hits; 00415 } 00416 //______________________________________________________________________________ 00417 TObjArray *StEventHelper::SelVertex(const char *sel,Int_t flag) 00418 { 00419 00420 TObjArray *conts = SelConts(sel); 00421 TObjArray *verts = new TObjArray(); 00422 Int_t ilast = conts->GetLast(); 00423 int nvtx =0; 00424 for (int idx=0;idx<=ilast;idx++) { 00425 StObjArray *arr = (StObjArray *)conts->At(idx); 00426 if (!arr) continue; 00427 int sz = arr->size(); 00428 if (!sz) continue; 00429 for (int ivx=0; ivx<sz; ivx++) { 00430 StVertex *vx = (StVertex*)arr->at(ivx); 00431 if (!vx) continue; 00432 if (vx->IsZombie()) continue; 00433 // See: StRoot/St_base/StObject.h also 00434 if ((flag&kMark2Draw) && !vx->TestBit(kMark2Draw)) continue; 00435 verts->Add(vx);nvtx++; 00436 } 00437 } 00438 delete conts; 00439 return verts; 00440 } 00441 //______________________________________________________________________________ 00442 TObjArray *StEventHelper::ExpandAndFilter(const TObject *eObj, int flag, TObjArray *out) 00443 { 00444 // eobj - TObject of StEvent objects (StHit,StTrack,StVertex or TObjArray of above) 00445 if (!out) {out = new TObjArray;} 00446 TObject *eobj = (TObject*)eObj; 00447 00448 int kind = Kind(eobj); 00449 if (kind&kHIT) {// input StHit 00450 if (!(flag&kHIT)) return out; 00451 int take=kind & (kUSE|kUNU|kFIT) &flag; 00452 if (take) out->Add(eobj); 00453 return out; 00454 }//endif StHit 00455 00456 //------------------------------------------------------------- 00457 if (kind&kTRK) {// input StTrack 00458 if (flag&kTRK) out->Add(eobj); 00459 if (!(flag&kHRR)) return out; 00460 StTrack *trk = (StTrack*)eobj; 00461 StTrackHelper trkh(trk); 00462 out->Add((TObject*)trkh.GetHits()); 00463 return out; 00464 }//endif StTrack 00465 00466 //------------------------------------------------------------- 00467 if (kind&kVTX) {// input StVertex 00468 if (flag&kVTX) out->Add(eobj); 00469 StVertex *vtx = (StVertex*)eobj; 00470 if (!(flag&(kTRK|kHRR))) return out; 00471 StVertexHelper vtxh(vtx); 00472 int n = vtxh.GetNTracks(); 00473 for (int i=-1;i<n;i++) { 00474 const TObject *to = vtxh.GetTrack(i); 00475 if (!to) continue; 00476 ExpandAndFilter(to,flag,out); 00477 } 00478 return out; 00479 }//endif StVertex 00480 //------------------------------------------------------------- 00481 00482 if (kind&kTRR) { // input TObjArray 00483 TObjArray *inp = (TObjArray *)eobj; 00484 inp->Compress(); 00485 int nbjs = inp->GetLast()+1; 00486 if (!nbjs) return 0; 00487 for (int i=0;i<nbjs;i++) { 00488 ExpandAndFilter(inp->At(i),flag,out); 00489 } 00490 return out; 00491 }// endif TObjArray 00492 00493 //------------------------------------------------------------- 00494 if (kind&kHRR) { // input HitsArray 00495 if (!(flag&kHRR)) return out; 00496 out->Add(eobj); 00497 return out; 00498 }// endif HitsArray 00499 return 0; 00500 } 00501 //______________________________________________________________________________ 00502 TObjArray *StEventHelper::MakePoints(TObjArray *inp, int flag) 00503 { 00504 static const Color_t plitra[]={kRed,kGreen,kBlue,kMagenta, kCyan}; 00505 static const int nlitra = sizeof(plitra)/sizeof(Color_t); 00506 int ilitra=0; 00507 inp->Compress(); 00508 int nbjs = inp->GetLast()+1; 00509 if (!nbjs) return 0; 00510 TObjArray *out = new TObjArray;out->SetOwner(); 00511 StPoints3DABC *p[3] ; int np=0; 00512 for (int i=0;i<nbjs;i++) { 00513 TObject *to = inp->At(i); 00514 int kind = Kind(to); 00515 if (!(kind&kHRR)) ilitra = (++ilitra)%nlitra; 00516 int take = (kind&flag); 00517 if (!take) continue; 00518 //?? if (take&kHIT) take -=kHIT; 00519 if (!take) continue; 00520 00521 np = 0; 00522 if (kind&kHIT) {np=1;p[0] = new StHitPoints ((StHit *)to );} 00523 else if (kind&kHRR && ((StPtrVecHit*)to)->size()) 00524 {np=1;p[0] = new StHitPoints ((StPtrVecHit*)to );} 00525 00526 else if (kind&kTRK) {np=3;p[0] = new StTrackPoints ((StTrack *)to ); 00527 p[1] = new StInnOutPoints((StTrack *)to,0); 00528 p[2] = new StInnOutPoints((StTrack *)to,1);} 00529 00530 else if (kind&kVTX) {np=1;p[0] = new StVertexPoints((StVertex *)to );} 00531 00532 for (int j=0;j<np;j++){p[j]->SetUniqueID(plitra[ilitra]); out->Add(p[j]);} 00533 }// 00534 00535 return out; 00536 } 00537 //______________________________________________________________________________ 00538 void StEventHelper::Break(int kase) 00539 { 00540 fprintf(stderr,"Break(%d)\n",kase); 00541 } 00542 00543 //______________________________________________________________________________ 00544 ClassImp(StPoints3DABC) 00545 void StPoints3DABC::Add(StPoints3DABC *add) 00546 { 00547 int n = add->fSize + fSize; 00548 if (n > fN) { 00549 if (n < fN*2) n = fN*2; 00550 Float_t *arr = new Float_t[n*3]; 00551 memcpy(arr, fXYZ, fSize*3*sizeof(Float_t)); 00552 delete [] fXYZ; fXYZ = arr; fN = n; 00553 } 00554 memcpy(fXYZ+fSize*3,add->fXYZ,add->fSize*3*sizeof(Float_t)); 00555 fSize+=add->fSize; 00556 } 00557 //______________________________________________________________________________ 00558 ClassImp(StTrackPoints) 00559 //______________________________________________________________________________ 00560 StTrackPoints::StTrackPoints(const StTrack *st,const char *name,const char *title) 00561 :StPoints3DABC(name,title,st) 00562 { 00563 Init(); 00564 } 00565 //______________________________________________________________________________ 00566 void StTrackPoints::Init() 00567 { 00568 if (fXYZ) return; 00569 StTrack *trk = ((StTrack*)fObj); 00570 StTrackHelper th(trk); 00571 fXYZ = th.GetPoints(fSize); 00572 fN = fSize; 00573 if (!fSize) { MakeZombie(); return;} 00574 } 00575 //______________________________________________________________________________ 00576 Int_t StTrackPoints::DistancetoPrimitive(Int_t px, Int_t py) 00577 { 00578 //*-*-*-*-*-*-*Compute distance from POINT px,py to a 3-D points *-*-*-*-*-*-* 00579 //*-* ===================================================== 00580 //*-* 00581 //*-* Compute the closest distance of approach from POINT px,py to each SEGMENT 00582 //*-* of the polyline. 00583 //*-* Returns when the distance found is below DistanceMaximum. 00584 //*-* The distance is computed in pixels units. 00585 //*-* 00586 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* 00587 00588 enum {inaxis = 7,mindist=10}; 00589 Float_t dist = 999999; 00590 00591 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin()); 00592 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin()); 00593 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax()); 00594 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax()); 00595 00596 TView *view = 0; 00597 //*-*- return if POINT is not in the user area 00598 if (px < puxmin - inaxis) goto END; 00599 if (py > puymin + inaxis) goto END; 00600 if (px > puxmax + inaxis) goto END; 00601 if (py < puymax - inaxis) goto END; 00602 00603 view = gPad->GetView(); 00604 if (!view) goto END; 00605 00606 {Int_t i; 00607 Float_t /*dpoint,*/alfa; 00608 Float_t xndc[3]; 00609 Int_t x1,y1,x0,y0; 00610 Int_t pointSize = fN*3; 00611 view->WCtoNDC(fXYZ, xndc); 00612 x0 = gPad->XtoAbsPixel(xndc[0]); 00613 y0 = gPad->YtoAbsPixel(xndc[1]); 00614 00615 float dif[2],difdif,cur[2],curcur,difcur; 00616 for (i=3;i<pointSize;i+=3) { 00617 view->WCtoNDC(fXYZ+i, xndc); 00618 x1 = gPad->XtoAbsPixel(xndc[0]); 00619 y1 = gPad->YtoAbsPixel(xndc[1]); 00620 dif[0] = x1-x0; dif[1]=y1-y0; 00621 cur[0] = x0-px; cur[1]=y0-py; 00622 difdif = (dif[0]*dif[0]+dif[1]*dif[1]); 00623 difcur = (dif[0]*cur[0]+dif[1]*cur[1]); 00624 curcur = cur[0]*cur[0]+cur[1]*cur[1]; 00625 if (difdif<mindist*mindist) { 00626 if ((i+3)<pointSize) continue; 00627 dist = curcur; break; 00628 } 00629 alfa = -difcur/difdif; 00630 00631 if (alfa<0.) {dist = curcur; break;} 00632 00633 x0=x1; y0=y1; 00634 if (alfa > 1.) { 00635 if (i+3 < pointSize) continue; 00636 dist = (px-x1)*(px-x1) + (py-y1)*(py-y1); break; 00637 } 00638 dist = curcur+alfa*(2*difcur+difdif*alfa); 00639 break; 00640 }} 00641 END: 00642 dist = TMath::Sqrt(dist); 00643 //VP if (dist <= mindist) { dist = 0; gPad->SetSelected(fObj);} 00644 if (dist <= mindist) { dist = 0; gPad->SetSelected(this);} 00645 00646 return Int_t(dist); 00647 } 00648 00649 00650 //______________________________________________________________________________ 00651 ClassImp(StVertexPoints) 00652 //______________________________________________________________________________ 00653 StVertexPoints::StVertexPoints(const StVertex *sv,const char *name,const char *title) 00654 :StPoints3DABC(name,title,sv) 00655 { 00656 SetBit(1); 00657 fSize = 1; fN =1; 00658 fXYZ = new Float_t[3]; 00659 fXYZ[0] = ((StVertex*)fObj)->position().x(); 00660 fXYZ[1] = ((StVertex*)fObj)->position().y(); 00661 fXYZ[2] = ((StVertex*)fObj)->position().z(); 00662 } 00663 //______________________________________________________________________________ 00664 ClassImp(StVertexPoints) 00665 //______________________________________________________________________________ 00666 StInnOutPoints::StInnOutPoints(const StTrack *st,int innout,const char *name,const char *title) 00667 :StPoints3DABC(name,title,st) 00668 { 00669 fSize = 1; fN =1; fInnOut=innout; 00670 const StTrackGeometry *geo = (fInnOut==0) ? st->geometry():st->outerGeometry(); 00671 fXYZ = new Float_t[3]; 00672 fXYZ[0] = geo->origin().x(); 00673 fXYZ[1] = geo->origin().y(); 00674 fXYZ[2] = geo->origin().z(); 00675 } 00676 00677 00678 ClassImp(StHitPoints) 00679 //______________________________________________________________________________ 00680 StHitPoints::StHitPoints(const StHit *sh,const char *name,const char *title) 00681 :StPoints3DABC(name,title,sh) 00682 { 00683 fSize = 1; fN =1; 00684 Init(); 00685 } 00686 //______________________________________________________________________________ 00687 StHitPoints::StHitPoints(const StRefArray *ar,const char *name,const char *title) 00688 :StPoints3DABC(name,title,ar) 00689 { 00690 fSize = ar->size(); 00691 fN = fSize; 00692 if (!fSize) return; 00693 fObj = (fSize==1) ? ar->front() : ar; 00694 Init(); 00695 } 00696 //______________________________________________________________________________ 00697 void StHitPoints::Init() 00698 { 00699 if (fXYZ) return; 00700 fXYZ = new Float_t[fN*3]; 00701 00702 int n=0; 00703 for (int i =0;i<fSize;i++) 00704 { 00705 StHit *hit= (fSize==1) ? (StHit*)fObj: (StHit*)((StRefArray*)fObj)->at(i); 00706 if (fSize>1 && !hit->trackReferenceCount()) continue; 00707 if (fSize>1 && !hit->usedInFit()) continue; 00708 StThreeVectorF v3 = hit->position(); 00709 fXYZ[n*3+0] = v3.x(); 00710 fXYZ[n*3+1] = v3.y(); 00711 fXYZ[n*3+2] = v3.z(); 00712 n++; 00713 } 00714 fN=n; fSize=n; 00715 } 00716 00717 //______________________________________________________________________________ 00718 ClassImp(StFilterABC) 00719 00720 int StFilterABC::fgDial=0; 00721 //______________________________________________________________________________ 00722 StFilterABC::StFilterABC(const char *name,bool active):TNamed(name,""),fActive(active) 00723 { 00724 } 00725 //______________________________________________________________________________ 00726 void StFilterABC::SetDefs() 00727 { 00728 for (int i=0;GetNams() && GetNams()[i]; i++) {GetPars()[i]=GetDefs()[i];} 00729 } 00730 00731 //______________________________________________________________________________ 00732 ClassImp(StFilterDef) 00733 StFilterDef::StFilterDef(const char *name,bool active):StFilterABC(name,active) 00734 { 00735 SetDefs(); 00736 00737 } 00738 //______________________________________________________________________________ 00739 const char **StFilterDef::GetNams() const 00740 { 00741 static const char *nams[] = { 00742 " RandomSelect ", 00743 " RxyMin ", 00744 " RxyMax ", 00745 " ZMin ", 00746 " ZMax ", 00747 " PhiMin ", 00748 " PhiMax ", 00749 " LenMin ", 00750 " LenMax ", 00751 " PtMin ", 00752 " PtMax ", 00753 " PseudoMin ", 00754 " PseudoMax ", 00755 " QMin ", 00756 " QMax ", 00757 " EncodedMethod", 00758 0}; 00759 return nams; 00760 } 00761 //______________________________________________________________________________ 00762 const float *StFilterDef::GetDefs() const 00763 { 00764 static const float defs[] = { 00765 /* RandomSelect=*/ 1.00, 00766 /* RxyMin =*/ 0.00, 00767 /* RxyMax =*/ 900.00, 00768 /* ZMin =*/ -900.00, 00769 /* ZMax =*/ +900.00, 00770 /* PhiMin =*/ -180.01, 00771 /* PhiMax =*/ +181.01, 00772 /* LenMin =*/ +0.00, 00773 /* LenMax =*/ +999.00, 00774 /* PtMin =*/ 0.00, 00775 /* PtMax =*/ 999.00, 00776 /* PseudoMin =*/ -999.00, 00777 /* PseudoMax =*/ 999.00, 00778 /* QMin =*/ -1 , 00779 /* QMax =*/ +1 , 00780 /* Encoded method*/ -1 , // The default value -1 menas all 00781 00782 0}; 00783 return defs; 00784 } 00785 00786 //______________________________________________________________________________ 00787 Int_t StFilterDef::Accept(StPoints3DABC *pnt,Color_t &color, Size_t&, Style_t&) 00788 { 00789 static TRandom rrr; 00790 float x,y,z,r2xy,phid,len,pt,ps,q; 00791 const TObject *to; 00792 const StTrack *trk; 00793 // set default color for tracks 00794 color = (((color-kRed)+1)%6)+kRed; 00795 00796 int cut = 1; 00797 if (fRandomSelect < 1. && fRandomSelect < rrr.Rndm())return 0; 00798 00799 z = pnt->GetZ(0); 00800 cut++; 00801 if (fZMin >z || z > fZMax) goto SKIP; 00802 00803 x = pnt->GetX(0); 00804 y = pnt->GetY(0); 00805 r2xy = x*x+y*y; 00806 cut++; 00807 if (fRxyMin*fRxyMin > r2xy || r2xy > fRxyMax*fRxyMax)goto SKIP; 00808 phid = atan2(y,x)*(180./M_PI); 00809 cut++; 00810 if (fPhiMin > phid || phid > fPhiMax) goto SKIP; 00811 to = pnt->GetObject(); 00812 if (!to) return 1; 00813 if (!to->InheritsFrom(StTrack::Class())) return 1; 00814 00815 // set default color for tracks 00816 // color = (((color-kRed)+1)%6)+kRed; 00817 00818 trk = (StTrack*)to; 00819 len = trk->length(); 00820 cut++; 00821 if (fLenMin >len || len > fLenMax) goto SKIP; 00822 pt = trk->geometry()->momentum().perp(); 00823 cut++; 00824 if (fPtMin >pt || pt > fPtMax) goto SKIP; 00825 ps = trk->geometry()->momentum().pseudoRapidity(); 00826 cut++; 00827 if (fPsMin >ps || ps > fPsMax) goto SKIP; 00828 q = trk->geometry()->charge(); 00829 cut++; 00830 if (fQMin >q || q > fQMax) goto SKIP; 00831 cut++; 00832 if ( (int(fEncodedMethod) != -1) && (trk->encodedMethod() != int(fEncodedMethod)) ) 00833 goto SKIP; 00834 return 1; 00835 00836 SKIP: return 0; 00837 00838 } 00839 00840 00841 //______________________________________________________________________________ 00842 ClassImp(StMuDstFilterHelper) 00843 StMuDstFilterHelper::StMuDstFilterHelper(const char *name,bool active):StFilterABC(name,active) 00844 { 00845 mBB = new BetheBloch(); 00846 SetDefs(); 00847 00848 } 00849 //______________________________________________________________________________ 00850 StMuDstFilterHelper::~StMuDstFilterHelper() 00851 { delete mBB;} 00852 //______________________________________________________________________________ 00853 const char **StMuDstFilterHelper::GetNams() const 00854 { 00855 static const char *nams[] = { 00856 " pCutHigh ", 00857 " nHitsCutHighP ", 00858 " pCutLow ", 00859 " nHitsCutLowP ", 00860 " chargeForLowP ", 00861 " dEdxMassCutHigh ", 00862 " dEdxFractionCutHigh ", 00863 " dEdxMassCutLow ", 00864 " dEdxFractionCutLow ", 00865 0 00866 }; 00867 return nams; 00868 } 00869 //______________________________________________________________________________ 00870 const float *StMuDstFilterHelper::GetDefs() const 00871 { 00872 static const float defs[] = { 00873 /* pCutHigh */ 2.0, // high momentum cut for RICH/Upsilon candidates 00874 /* nHitsCutHighP */ 10, // nHits cut for all tracks 00875 /* pCutLow */ 0.2, // low momentum cut 00876 /* nHitsCutLowP */ 15, 00877 /* chargeForLowP */ -1, // charge for tracks with pCutLow < p < pCutHigh, set to 0 for all tracks 00878 /* dEdxMassCutHigh */ 0.939, // cut below BetheBloch(p/dEdxMassCutHigh), e.g. proton-band 00879 /* dEdxFractionCutHigh */ 0.6, // cut fraction of dEdx-band, i.e. dEdxFractionCut * BetheBloch(p/dEdxMassCut) 00880 /* dEdxMassCutLow */ 0.494, // cut above BetheBloch(p/dEdxMassCutLow), e.g. kaon-band 00881 /* dEdxFractionCutLow */ 1.1, 00882 0 00883 }; 00884 return defs; 00885 } 00886 //______________________________________________________________________________ 00887 Int_t StMuDstFilterHelper::Accept(const StTrack* track) { 00888 00889 float pCutHigh = fpCutHigh; // high momentum cut for RICH/Upsilon candidates 00890 int nHitsCutHighP = int(fnHitsCutHighP); // nHits cut for all tracks 00891 00892 // following cuts apply only for tracks with pCutLow < p <pHigh 00893 float pCutLow = fpCutLow; // low momentum cut 00894 int nHitsCutLowP = int(fnHitsCutLowP); 00895 int chargeForLowP = int(fchargeForLowP); // charge for tracks with pCutLow < p < fpCutHigh, set to 0 for all tracks 00896 float dEdxMassCutHigh = fdEdxMassCutHigh; // cut below BetheBloch(p/dEdxMassCutHigh), e.g. proton-band 00897 float dEdxFractionCutHigh = fdEdxFractionCutHigh;// cut fraction of dEdx-band, i.e. dEdxFractionCut * BetheBloch(p/dEdxMassCut) 00898 float dEdxMassCutLow = fdEdxMassCutLow; // cut above BetheBloch(p/dEdxMassCutLow), e.g. kaon-band 00899 float dEdxFractionCutLow = fdEdxFractionCutLow; 00900 00901 int iret = 0; 00902 int chargeOK = 0; 00903 int dedxOK = 0; 00904 00905 float magnitude = track->geometry()->momentum().magnitude(); 00906 int nPoints = track->detectorInfo()->numberOfPoints(); 00907 00908 if ( magnitude > pCutHigh && nPoints >= nHitsCutHighP) iret = 1; 00909 else { 00910 if ( magnitude > pCutLow && nPoints >= nHitsCutLowP ) 00911 { 00912 // check charge 00913 if (chargeForLowP==0) 00914 chargeOK = 1; 00915 else if (track->geometry()->charge() == chargeForLowP) 00916 chargeOK = 1; 00917 00918 // check dEdx 00919 // if (mBB==0) mBB = new BetheBloch(); 00920 float dedxHigh = dEdxFractionCutHigh * mBB->Sirrf(magnitude/dEdxMassCutHigh); 00921 float dedxLow = dEdxFractionCutLow * mBB->Sirrf(magnitude/dEdxMassCutLow); 00922 float dedx = 0; 00923 00924 // get track dEdx 00925 const StSPtrVecTrackPidTraits& traits = track->pidTraits(); 00926 StDedxPidTraits* dedxPidTr; 00927 for (unsigned int itrait = 0; itrait < traits.size(); itrait++){ 00928 dedxPidTr = 0; 00929 if (traits[itrait]->detector() == kTpcId) { 00930 StTrackPidTraits* thisTrait = traits[itrait]; 00931 dedxPidTr = dynamic_cast<StDedxPidTraits*>(thisTrait); 00932 if (dedxPidTr && dedxPidTr->method() == kTruncatedMeanId) { 00933 // adjust L3 dE/dx by a factor of 2 to match offline 00934 dedx = 2 * dedxPidTr->mean(); 00935 } 00936 } 00937 } 00938 if (dedx > dedxHigh && dedx > dedxLow) 00939 dedxOK = 1; 00940 // final answer 00941 iret = chargeOK * dedxOK; 00942 } // if (pCutLow && nHitsCutLowP) 00943 } 00944 return iret; 00945 } 00946 00947 //______________________________________________________________________________ 00948 Int_t StMuDstFilterHelper::Accept(StPoints3DABC *pnt) 00949 { 00950 const TObject *to; 00951 const StTrack *trk; 00952 to = pnt->GetObject(); 00953 if (!to) return 1; 00954 if (!to->InheritsFrom(StTrack::Class())) return 1; 00955 trk = (StTrack*)to; 00956 return Accept(trk); 00957 } 00958 00959 //______________________________________________________________________________ 00960 ClassImp(StColorFilterHelper) 00961 StColorFilterHelper::StColorFilterHelper(const char *name,bool active):StFilterABC(name,active) 00962 { 00963 fPidAlgorithm = new StTpcDedxPidAlgorithm();; 00964 fElectron = StElectron::instance(); 00965 fPion = StPionPlus::instance(); 00966 fKaon = StKaonPlus::instance(); 00967 fProton = StProton::instance(); 00968 00969 SetDefs(); 00970 00971 } 00972 //______________________________________________________________________________ 00973 StColorFilterHelper::~StColorFilterHelper() 00974 { delete fPidAlgorithm;} 00975 //______________________________________________________________________________ 00976 const char **StColorFilterHelper::GetNams() const 00977 { 00978 static const char *nams[] = { 00979 " Electron sigma ", 00980 " Electron color ", 00981 " Pion sigma ", 00982 " Pion color ", 00983 " Kaon sigma ", 00984 " Kaon color ", 00985 " Proton sigma ", 00986 " Proton color ", 00987 " others sigma ", 00988 " others color ", 00989 0 00990 }; 00991 return nams; 00992 } 00993 //______________________________________________________________________________ 00994 const float *StColorFilterHelper::GetDefs() const 00995 { 00996 static const float defs[] = { 00997 /* fNSigmaElectron*/ 1 , // nSigma cut for electron 00998 /* fNColorElectron*/ 2 , // the color index for electron 00999 /* fNSigmaPion */ 1 , // nSigma cut for electron 01000 /* fNColorPion */ 3 , // the color index for pion 01001 /* fNSigmaKaon */ 1 , // nSigma cut for pion 01002 /* fNColorKaon */ 4 , // the color index for kaon 01003 /* fNSigmaProton */ 1 , // nSigma cut for kaon 01004 /* fNColorProton */ 6 , // the color index kaon 01005 01006 /* fNSigmaOther*/ -1, // nSigma cut for other types 01007 /* fNColorOther*/ 0, // the color index for other types 01008 0 01009 }; 01010 return defs; 01011 } 01012 //______________________________________________________________________________ 01013 Int_t StColorFilterHelper::Accept(const StTrack* track, Color_t &color, Size_t&size, Style_t&) { 01014 01015 float sigmaElectron = fNSigmaElectron ; // nSigna cut for electron 01016 Color_t colorElectron = (Color_t)fNColorElectron ; // the color index for electron 01017 01018 float sigmaPion = fNSigmaPion ; // nSigna cut for electron 01019 Color_t colorPion = (Color_t)fNColorPion ; // the color index for pion 01020 01021 float sigmaKaon = fNSigmaKaon ; // nSigna cut for pion 01022 Color_t colorKaon = (Color_t)fNColorKaon ; // the color index for kaon 01023 01024 float sigmaProton = fNSigmaProton ; // nSigna cut for kaon 01025 Color_t colorProton = (Color_t)fNColorProton ; // the color index kaon 01026 01027 // float sigmaOther = fNSigmaOther ; // nSigna cut for other types 01028 Color_t colorOther = (Color_t)fNColorOther ; // the color index for other types 01029 01030 // Fisyak's color schema 01031 01032 /* const StParticleDefinition* pd = */ track->pidTraits(*fPidAlgorithm); 01033 01034 color = colorOther; 01035 size = 1; 01036 01037 if (TMath::Abs(fPidAlgorithm->numberOfSigma(fElectron)) < sigmaElectron) 01038 { color = colorElectron; size = 2; } 01039 01040 if (TMath::Abs(fPidAlgorithm->numberOfSigma(fKaon)) < sigmaKaon) 01041 { color = colorKaon; size = 4; } 01042 01043 if (TMath::Abs(fPidAlgorithm->numberOfSigma(fPion)) < sigmaPion) 01044 { color = colorPion; size = 5; } 01045 01046 if (TMath::Abs(fPidAlgorithm->numberOfSigma(fProton)) < sigmaProton) 01047 { color = colorProton; size = 3; } 01048 01049 return 1; 01050 } 01051 01052 //______________________________________________________________________________ 01053 Int_t StColorFilterHelper::Accept(StPoints3DABC *pnt, Color_t&color, Size_t&size, Style_t&style) 01054 { 01055 const TObject *to; 01056 const StTrack *trk; 01057 to = pnt->GetObject(); 01058 if (!to) return 1; 01059 if (!to->InheritsFrom(StTrack::Class())) return 1; 01060 trk = (StTrack*)to; 01061 return Accept(trk,color,size,style); 01062 } 01063 //______________________________________________________________________________ 01064 //______________________________________________________________________________ 01065 //______________________________________________________________________________ 01066 ClassImp(StVertexHelper) 01067 StVertexHelper::StVertexHelper(const StVertex *vtx) 01068 { SetVertex(vtx);} 01069 StVertexHelper::StVertexHelper(const StEvent *evt) 01070 { SetVertex(evt->primaryVertex(0));} 01071 //______________________________________________________________________________ 01072 void StVertexHelper::SetVertex(const StVertex *vtx){fVtx = vtx;} 01073 int StVertexHelper::GetType() {return (int)fVtx->type();} 01074 int StVertexHelper::GetFlag() {return fVtx->flag();}; 01075 int StVertexHelper::GetNTracks() {return fVtx->numberOfDaughters();} 01076 //______________________________________________________________________________ 01077 const StThreeVectorF &StVertexHelper::GetPoint() 01078 { 01079 return fVtx->position(); 01080 } 01081 //______________________________________________________________________________ 01082 const StTrack *StVertexHelper::GetTrack(int idx) // -1=parent track 01083 { 01084 if (idx==-1) return fVtx->parent(); 01085 if (idx>= GetNTracks()) return 0; 01086 return fVtx->daughter((UInt_t)idx); 01087 } 01088 //______________________________________________________________________________ 01089 const float *StVertexHelper::GetErrMtx() 01090 { 01092 01093 StMatrixF mxF = fVtx->covariantMatrix(); 01094 int jj=0; 01095 for (int i=0;i< 3;i++) { 01096 for (int j=0;j<=i;j++) { 01097 fErrMtx[jj++] = mxF(i+1,j+1);}} 01098 return fErrMtx; 01099 } 01100 //______________________________________________________________________________ 01101 //______________________________________________________________________________ 01102 //______________________________________________________________________________ 01103 ClassImp(StTrackHelper) 01104 StTrackHelper::StTrackHelper(const StTrack *trk) 01105 { 01106 fHelx[0]=0; fHelx[1]=0; 01107 fTHlx[0]=0; fTHlx[1]=0; 01108 SetTrack(trk); 01109 } 01110 01111 StTrackHelper::~StTrackHelper() 01112 { 01113 delete fHelx[0];delete fHelx[1]; 01114 delete fTHlx[0];delete fTHlx[1]; 01115 } 01116 void StTrackHelper::SetTrack(const StTrack *trk){fTrk=trk;fHits=0;GetNHits();} 01117 int StTrackHelper::GetType() const {return fTrk->type();} 01118 int StTrackHelper::GetFlag() const {return fTrk->flag();} 01119 int StTrackHelper::GetCharge() const {return fTrk->geometry()->charge();} 01120 const StVertex *StTrackHelper::GetParent() const {return fTrk->vertex();} 01121 float StTrackHelper::GetImpact() const {return fTrk->impactParameter();} 01122 float StTrackHelper::GetCurv() const {return GetTHelix(0)->GetRho() ;} 01123 float StTrackHelper::GetLength() const {return fTrk->length();} 01124 const StThreeVectorF &StTrackHelper::GetFirstPoint() const {return fTrk->geometry()->origin();} 01125 const StThreeVectorF &StTrackHelper::GetLastPoint() const {return fTrk->outerGeometry()->origin();} 01126 const StThreeVectorF &StTrackHelper::GetMom() const {return fTrk->geometry()->momentum();} 01127 //______________________________________________________________________________ 01128 StPhysicalHelixD *StTrackHelper::GetHelix(int idx) const 01129 { 01130 if (!fHelx[idx]) fHelx[idx]= new StPhysicalHelixD; 01131 *fHelx[idx] = (idx==0) ? fTrk->geometry()->helix():fTrk->outerGeometry()->helix(); 01132 return fHelx[idx]; 01133 } 01134 //______________________________________________________________________________ 01135 THelixTrack *StTrackHelper::GetTHelix(int idx) const 01136 { 01137 StPhysicalHelixD *hlx = GetHelix(idx); 01138 fTHlx[idx] = StEventHelper::MyHelix(fTHlx[idx],hlx); 01139 return fTHlx[idx]; 01140 } 01141 //______________________________________________________________________________ 01142 const StPtrVecHit *StTrackHelper::GetHits() const 01143 { 01144 if (fHits) return fHits; 01145 const StTrackDetectorInfo *tdi = fTrk->detectorInfo(); 01146 if (!tdi) return 0; 01147 fHits = &tdi->hits(); 01148 return fHits; 01149 } 01150 //______________________________________________________________________________ 01151 int StTrackHelper::GetNHits() const 01152 { 01153 if (fHits) return fHits->size(); 01154 GetHits(); 01155 return (fHits)? fHits->size():0; 01156 } 01157 //______________________________________________________________________________ 01158 const StHit *StTrackHelper::GetHit(int idx) const 01159 { 01160 if (idx<0) return 0; 01161 if (idx>=GetNHits()) return 0; 01162 if (!fHits) return 0; 01163 return fHits->at(idx); 01164 } 01165 //______________________________________________________________________________ 01166 int StTrackHelper::numberOfFitPoints(int det) const 01167 { 01168 const StTrackFitTraits& trait = fTrk->fitTraits(); 01169 return (det)? trait.numberOfFitPoints((StDetectorId)det): trait.numberOfFitPoints(); 01170 } 01171 //______________________________________________________________________________ 01172 StMCTruth StTrackHelper::GetTruth(int byNumb,double rXYMin,double rXYMax) const 01173 { 01174 StMCPivotTruth pivo(1); 01175 int nHits = GetNHits(); 01176 int nUsed=0; 01177 for (int jh=0;jh<nHits;jh++) { 01178 const StHit *hit = GetHit(jh); 01179 double r = sqrt(pow(hit->position().x(),2)+pow(hit->position().y(),2)); 01180 if (r<rXYMin) continue; 01181 if (r>rXYMax) continue; 01182 int idTruth=hit->idTruth(); 01183 int wtTruth=hit->qaTruth(); 01184 if (!wtTruth) wtTruth=1; 01185 // if (!idTruth || !wtTruth) { 01186 // Warning("GetTruth","idTruth,wtTruth= %d %d",idTruth,wtTruth); 01187 // continue;} 01188 nUsed++; pivo.Add(idTruth,wtTruth); 01189 } 01190 if (!nUsed) return 0; 01191 return pivo.Get(byNumb); 01192 } 01193 01194 //______________________________________________________________________________ 01195 Float_t *StTrackHelper::GetPoints(int &npoints) const 01196 { 01197 static int ndebug=0; ndebug++; 01198 npoints=0; 01199 double len,len0,len1; 01200 len = fTrk->length(); 01201 if (len <= 0.0001) { 01202 Warning("GetPoints","Zero length %s(%p), IGNORED",fTrk->ClassName(),(void*)fTrk); 01203 return 0; 01204 } 01205 01206 GetHelix(0); GetHelix(1); 01207 for (int i=0;i<2;i++) {fTHlx[i] = StEventHelper::MyHelix(fTHlx[i],fHelx[i]);} 01208 01209 len0 = fTHlx[0]->Step(fTHlx[1]->GetXYZ()); 01210 double rho0 = fTHlx[0]->GetRho(); 01211 double rho1 = fTHlx[1]->GetRho(); 01212 double drho = (rho1-rho0)/(len0*fTHlx[0]->GetCos()); 01213 fTHlx[0]->Set(rho0,drho); 01214 fTHlx[1]->Set(rho1,drho); 01215 fTHlx[1]->Backward(); 01216 npoints = abs(int(len*fTHlx[0]->GetCos()*rho0*90))+2; 01217 double step = 1./(npoints-1); 01218 len0 = fTHlx[0]->Step(fTHlx[1]->GetXYZ()); 01219 len1 = fTHlx[1]->Step(fTHlx[0]->GetXYZ()); 01220 float *arr = new Float_t[npoints*3]; 01221 double xyz[3][3]; 01222 for (int i =0;i<npoints;i++) 01223 { 01224 double s0 = i*step; 01225 double s1 = 1.-s0; 01226 fTHlx[0]->Step(s0*len0,xyz[0]); 01227 fTHlx[1]->Step(s1*len1,xyz[1]); 01228 s0 = s0*s0*s0; s1 = s1*s1*s1; 01229 double tmp = s0+s1; 01230 s0 /=tmp; s1 /=tmp; 01231 TCL::vlinco(xyz[0],s1,xyz[1],s0,xyz[2],3); 01232 TCL::ucopy(xyz[2],arr+i*3,3); 01233 } 01234 return arr; 01235 } 01236 01237 01238 //______________________________________________________________________________ 01239 //______________________________________________________________________________ 01240 //______________________________________________________________________________ 01241 ClassImp(StHitHelper) 01242 StHitHelper::StHitHelper(const StHit *hit){fHit = hit;} 01243 void StHitHelper::SetHit(const StHit *hit) {fHit = hit;} 01244 int StHitHelper::GetDetId() {return fHit->detector();} 01245 int StHitHelper::GetFlag() {return fHit->flag();} 01246 float StHitHelper::GetCharge() {return fHit->charge();} 01247 int StHitHelper::IsUsed() {return fHit->trackReferenceCount();} 01248 int StHitHelper::IsFit() {return fHit->usedInFit();} 01249 const StThreeVectorF &StHitHelper::GetPoint() {return fHit->position();} 01250 //______________________________________________________________________________ 01251 //______________________________________________________________________________ 01252 //______________________________________________________________________________ 01253 01254 ClassImp(StErrorHelper) 01255 //_____________________________________________________________________________ 01256 StErrorHelper::StErrorHelper() 01257 { 01258 fNErr=0; fNTot=0; fKErr=0; 01259 fMap = new TExMap; 01260 fArr = new TArrayI; 01261 } 01262 //_____________________________________________________________________________ 01263 StErrorHelper::~StErrorHelper() 01264 { 01265 delete fMap; 01266 delete fArr; 01267 } 01268 //_____________________________________________________________________________ 01269 void StErrorHelper::Add(int errn) 01270 { 01271 fNTot++; 01272 if(!errn) return; 01273 fNErr++; 01274 (*fMap)(errn)++; 01275 } 01276 //_____________________________________________________________________________ 01277 void StErrorHelper::MakeArray() 01278 { 01279 if (!fNErr) return; 01280 fKErr = fMap->GetSize(); 01281 fArr->Set(fKErr*3); 01282 TExMapIter it(fMap); 01283 long lerr,lnum; 01284 int idx=0; 01285 while(it.Next(lerr,lnum)) { 01286 (*fArr)[idx+ 0] = lnum; 01287 (*fArr)[idx+fKErr] = lerr; 01288 idx++; 01289 } 01290 01291 TMath::Sort(fKErr, fArr->GetArray(), fArr->GetArray()+2*fKErr); 01292 } 01293 //_____________________________________________________________________________ 01294 void StErrorHelper::Print(const char* txt) const 01295 { 01296 StErrorHelper *This = (StErrorHelper *)this; 01297 This->MakeArray(); 01298 if (!txt) txt=""; 01299 printf("StEvent Error Summary:%s\n",txt); 01300 01301 printf("%4d -%8d(%4d)\n",0,0,fNTot-fNErr); 01302 int *nrr=fArr->GetArray(); 01303 int *krr=nrr+fKErr; 01304 int *idx=krr+fKErr; 01305 for (int i=0;i<fKErr;i++) { 01306 int j = idx[i]; 01307 printf("%4d -%8d(%4d) //%s\n",i+1,krr[j],nrr[j],Say(krr[j]).Data()); 01308 } 01309 } 01310 //_____________________________________________________________________________ 01311 01312 TString StErrorHelper::Say(int ierr,const char *klass) 01313 { 01314 static const char *TabErr[] = 01315 { 01316 "StTrack" ,"mFlag" ,"1","2","is Negative", 01317 "StTrack" ,"mFlag" ,"1","3","is Zero", 01318 01319 "StTrack" ,"mImpactParameter" ,"2","1","is NaN", 01320 "StTrack" ,"mImpactParameter" ,"2","2","is huge", 01321 01322 "StTrack" ,"mLength" ,"3","1","is NaN", 01323 "StTrack" ,"mLength" ,"3","2","is huge", 01324 "StTrack" ,"mLength" ,"3","3","is too small", 01325 "StTrack" ,"mLength" ,"3","4","contradicts to In/Out distance", 01326 "StTrack" ,"mLength" ,"3","5","helix out of Zmax", 01327 "StTrack" ,"mLength" ,"3","6","helix out of Rmax", 01328 01329 "StTrack" ,"mGeometry" ,"4","2","iz zero", 01330 "StTrack" ,"mGeometry" ,"4","0","StTrackGeometry", 01331 01332 "StTrack" ,"mOuterGeometry" ,"5","2","iz zero", 01333 "StTrack" ,"mOuterGeometry" ,"5","0","StTrackGeometry", 01334 01335 "StTrack" ,"mDetectorInfo" ,"6","2","iz zero", 01336 "StTrack" ,"mDetectorInfo" ,"6","0","StTrackDetectorInfo", 01337 01338 "StTrackGeometry" ,"Helix" ,"1","0","StPhysicalHelixD", 01339 "StTrackGeometry" ,"Helix" ,"1","2","out of zMax", 01340 "StTrackGeometry" ,"Helix" ,"1","3","out of rMax", 01341 01342 "StTrackDetectorInfo" ,"mFirstPoint" ,"1","0","StThreeVectorF", 01343 "StTrackDetectorInfo" ,"mFirstPoint" ,"1","2","out of zMax", 01344 "StTrackDetectorInfo" ,"mFirstPoint" ,"1","3","out of rMax", 01345 01346 "StTrackDetectorInfo" ,"mLastPoint" ,"2","0","StThreeVectorF", 01347 "StTrackDetectorInfo" ,"mLastPoint" ,"2","2","out of zMax", 01348 "StTrackDetectorInfo" ,"mFLastPoint" ,"2","3","out of rMax", 01349 01350 01351 "StPhysicalHelixD" ,"mDipAngle" ,"1","1","is NaN", 01352 "StPhysicalHelixD" ,"mDipAngle" ,"1","2","> Py/2", 01353 "StPhysicalHelixD" ,"mDipAngle" ,"1","3","== Py/2", 01354 01355 "StPhysicalHelixD" ,"mCurvature" ,"2","1","is NaN", 01356 "StPhysicalHelixD" ,"mCurvature" ,"2","2","too big", 01357 "StPhysicalHelixD" ,"mCurvature" ,"2","3","is Negaive", 01358 01359 "StPhysicalHelixD" ,"mOrigin" ,"3","0","StThreeVectorD", 01360 "StPhysicalHelixD" ,"mH" ,"4","2","!= 1 or -1", 01361 01362 "StThreeVectorD" ,"mX1" ,"1","1","is NaN", 01363 "StThreeVectorD" ,"mX1" ,"1","2","too big", 01364 "StThreeVectorD" ,"mX2" ,"2","1","is NaN", 01365 "StThreeVectorD" ,"mX2" ,"2","2","too big", 01366 "StThreeVectorD" ,"mX3" ,"3","1","is NaN", 01367 "StThreeVectorD" ,"mX3" ,"3","2","too big", 01368 01369 01370 "StThreeVectorF" ,"mX1" ,"1","1","is NaN", 01371 "StThreeVectorF" ,"mX1" ,"1","2","too big", 01372 "StThreeVectorF" ,"mX2" ,"2","1","is NaN", 01373 "StThreeVectorF" ,"mX2" ,"2","2","too big", 01374 "StThreeVectorF" ,"mX3" ,"3","1","is NaN", 01375 "StThreeVectorF" ,"mX3" ,"3","2","too big", 01376 0}; 01377 TString ts; 01378 01379 int jmm = ierr%10; 01380 int jrr = (ierr/10)%10; 01381 for (const char **jt=TabErr;*jt;jt+=5) { 01382 if (strcmp(klass,*jt) ) continue; 01383 if (atoi(jt[2]) != jmm) continue; 01384 if (atoi(jt[3]) != jrr) continue; 01385 ts+=jt[1]; 01386 if (jrr) { ts+=": "; ts+=jt[4];} 01387 else { ts+="."; ts+=Say(ierr/100,jt[4]);} 01388 return ts; 01389 } 01390 ts="***Unknown***"; 01391 return ts; 01392 } 01393 //________________________________________________________________________________ 01394 //________________________________________________________________________________ 01395 StHitIter::StHitIter() 01396 { 01397 fEvent = 0; 01398 Reset(); 01399 } 01400 //________________________________________________________________________________ 01401 void StHitIter::Reset() 01402 { 01403 memset(fO, 0,sizeof(fO)); 01404 memset(fN, 0,sizeof(fN)); 01405 memset(fI ,-1,sizeof(fI)); 01406 } 01407 //________________________________________________________________________________ 01408 StHitIter &StHitIter::Next() 01409 { 01410 int kase = 0; 01411 fO[0]= 0; 01412 01413 while(2005) { 01414 if (++fI[kase] >= fN[kase]) {kase++;continue;} 01415 fO[kase] = GetO(kase+1,fI[kase]); 01416 if (!fO[kase]) break; 01417 if (!kase ) break; 01418 fN[kase-1] = GetN(kase); 01419 fI[kase-1] = -1; 01420 kase--; 01421 } 01422 return *this; 01423 } 01424 //________________________________________________________________________________ 01425 //________________________________________________________________________________ 01426 //________________________________________________________________________________ 01427 StSvtHitIter::StSvtHitIter(StEvent *ev):StHitIter() 01428 { 01429 fEvent = ev; 01430 fO[4] = fEvent->svtHitCollection(); 01431 fN[3] = 1; 01432 ++(*this); 01433 } 01434 //________________________________________________________________________________ 01435 StObject *StSvtHitIter::GetO(int lev,int idx) 01436 { 01437 switch (lev) { 01438 01439 case 1: //StSvtHit 01440 return ((StSvtWaferHitCollection *)fO[1])->hits()[idx]; 01441 case 2: //Wafer 01442 return ((StSvtLadderHitCollection*)fO[2])->wafer (idx); 01443 case 3: //Ladder 01444 return ((StSvtBarrelHitCollection*)fO[3])->ladder(idx); 01445 case 4: //Ladder 01446 return ((StSvtHitCollection *)fO[4])->barrel(idx); 01447 case 5: //SvtHitCollection 01448 return fEvent->svtHitCollection(); 01449 default: return 0; 01450 } 01451 } 01452 01453 //________________________________________________________________________________ 01454 int StSvtHitIter::GetN(int lev) const 01455 { 01456 switch (lev) { 01457 case 1: //N StSvtHits 01458 return ((StSvtWaferHitCollection *)fO[1])->hits().size(); 01459 case 2: //N Wafers 01460 return ((StSvtLadderHitCollection*)fO[2])->numberOfWafers(); 01461 case 3: //N Ladders 01462 return ((StSvtBarrelHitCollection*)fO[3])->numberOfLadders(); 01463 case 4: //N Barrels 01464 return ((StSvtHitCollection *)fO[4])->numberOfBarrels(); 01465 case 5: //N SvtHitCollections 01466 return 1; 01467 default: return 0; 01468 } 01469 } 01470 01471 01472

Generated on Sun Mar 15 04:54:28 2009 for StRoot by doxygen 1.3.7