StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
EventT.cxx
1 // ROOT headers
2 #include "TH1.h"
3 #include "TH2.h"
4 #include "TStyle.h"
5 #include "TCanvas.h"
6 #include "TF1.h"
7 #include "TGeoMatrix.h"
8 #include "TKey.h"
9 #include "TDirectory.h"
10 #include "TClass.h"
11 #include "TRVector.h"
12 #include "TRSymMatrix.h"
13 
14 // StRoot headers
15 #include "StGmtHit.h"
16 #include "StGmtHitCollection.h"
17 #include "StGmtCollection.h"
18 #include "StarRoot/THelixTrack.h"
19 #include "EventT.h"
20 #include "TrackT.h"
21 #include "HitT.h"
22 #include "StEvent.h"
23 #include "StPrimaryVertex.h"
24 #include "StEventInfo.h"
25 #include "StEventSummary.h"
26 #include "StTrack.h"
27 #include "StTrackNode.h"
28 #include "StPrimaryTrack.h"
29 #include "StGlobalTrack.h"
30 #include "StTrackDetectorInfo.h"
31 #include "StTrackGeometry.h"
32 #include "StDcaGeometry.h"
33 #include "St_base/StMessMgr.h"
34 
35 // C++ headers
36 #include <cassert>
37 
38 #define __OnlyGlobal__
39 
40 TClonesArray *EventT::fgTracks = 0;
41 TClonesArray *EventT::fgHits = 0;
42 THashList *EventT::fRotList = 0;
43 
44 static int _debug = 0;
45 
46 //_________________
47 EventT::EventT() : fIsValid(kFALSE) {
48  // Create an EventT object.
49  // When the constructor is invoked for the first time, the class static
50  // variable fgTracks is 0 and the TClonesArray fgTracks is created.
51 
52  if (!fgTracks) fgTracks = new TClonesArray("TrackT", 1000);
53  fTracks = fgTracks;
54  fNtrack = 0;
55  if (!fgHits) fgHits = new TClonesArray("HitT", 1000);
56  fHits = fgHits;
57  fNhit = 0;
58 }
59 
60 //_________________
61 EventT::~EventT() {
62  Clear();
63  SafeDelete(fRotList);
64 }
65 
66 //_________________
67 Int_t EventT::Build(StEvent *pEventT, Double_t pCut) {
68 
69 #ifdef Extra
70  static const Int_t NoFitPointCutForGoodTrackT = 15;
71 #endif
72  Clear();
73  Int_t iok = 1;
74  fIsValid = kFALSE;
75  if (! pEventT) return iok;
76  StGmtCollection* GmtCollection = pEventT->gmtCollection();
77  if (! GmtCollection) {
78  LOG_INFO << "No GMT Collections" << endm;
79  return iok;
80  }
81  StSPtrVecTrackNode& theNodes = pEventT->trackNodes();
82  UInt_t nnodes = theNodes.size();
83  if (! nnodes) {
84  LOG_INFO << "No tracks" << endm;
85  return iok;
86  }
87  StEventInfo* info = pEventT->info();
88  Int_t ev = 0, run = 0, time = 0;
89  if (info) {
90  ev = info->id();
91  run = info->runId();
92  time = info->time();
93  }
94  StEventSummary* summary = pEventT->summary();
95  Double32_t field = 0;
96  if (summary) field = summary->magneticField();
97  SetHeader(ev,run,time,field);
98  SetFlag(1);
99  Int_t nTotMatch = 0;
100  // Create and Fill the TrackT objects
101  for (UInt_t i=0; i<nnodes; i++) {
102  StTrackNode *node = theNodes[i];
103  if (! node) continue;
104  StTrack *Track = node->track(global);
105  if (! Track) continue;
106  StGlobalTrack *gTrack = (StGlobalTrack *) Track;
107  StDcaGeometry *dca = gTrack->dcaGeometry();
108  if (! dca) continue;
109  StPrimaryTrack *pTrack = (StPrimaryTrack *) node->track(primary);
110 #ifndef __OnlyGlobal__
111  if (! pTrack) continue;
112  Track = (StTrack *) pTrack;
113 #endif
114  StThreeVectorD g3 = Track->geometry()->momentum();
115  Double_t pT = g3.perp();
116  Double_t pMom = g3.mag();
117  if (pMom < pCut) continue;
118  TrackT *track = AddTrackT();
119  Double_t InvpT = 0;
120  Double_t TanL = 999999;
121  if (TMath::Abs(pT) > 1.e-7) {
122  InvpT = Track->geometry()->charge()/pT;
123  TanL = g3.z()/pT;
124  }
125  track->SetInvpT(InvpT);
126  track->SetPhi(TMath::ATan2(g3.y(),g3.x()));
127  track->SetTanL(TanL);
128  static const Double_t EC = 2.9979251E-4;
129  Double_t Rho = - EC*InvpT*field;
130  track->SetRho(Rho);
131  track->SetLength(Track->length());
132  StTrackDetectorInfo* dinfo=Track->detectorInfo();
133  track->SetNpoint(dinfo->numberOfPoints());
134  track->SetNPpoint(Track->numberOfPossiblePoints());
135  track->SetN(0);
136  StPhysicalHelixD helixO = Track->outerGeometry()->helix();
137  Double_t R_tof[2]= {210., 216.}; // inner and outer surfaces
138  pair<Double_t,Double_t> shR = helixO.pathLength(0.5*(R_tof[0]+R_tof[1]));
139 
140  if (TMath::Abs(shR.first) > 200 && TMath::Abs(shR.second) > 200) continue;
141 
142  Double_t stepR = (shR.first > 0) ? shR.first : shR.second;
143  StThreeVectorD xyzR = helixO.at(stepR);
144  Double_t phiR = TMath::RadToDeg()*xyzR.phi();
145  if (_debug) {
146  LOG_INFO << "\t shR " << shR.first << "\t" << shR.second << "\tstepR " << stepR
147  << "\txyzR\t" << xyzR << "\tphiR\t" << phiR << endm;
148  }
149  Int_t NoHitPerTrack = 0;
150  THashList *fRotList = RotMatrices();
151  if (! fRotList) continue;
152  TIter next(fRotList);
153  TGeoHMatrix *comb = 0;
154 
155  while ((comb = (TGeoHMatrix *) next())) {
156  TString combName(comb->GetName());
157  if (! combName.BeginsWith("R")) continue;
158  Int_t Id;
159  sscanf(comb->GetName()+1,"%i",&Id);
160  UInt_t module = Id;
161  StGmtHitCollection* GmtHitCollection = GmtCollection->getHitCollection(module);
162  if (! GmtHitCollection) { cout << "No GMT HitT Collections for mudule " << module << endl; continue;}
163  StSPtrVecGmtHit& hitvec = GmtHitCollection->getHitVec();
164  UInt_t NoHits = hitvec.size();
165  if (! NoHits) continue;
166  if (_debug) {
167  LOG_INFO << comb->GetName() << "\tmodule = " << module << endm;
168  comb->Print();
169  }
170  static Double_t dz[2] = {50.00, 2.10};
171  static Double_t dx[2] = {50.00, 3.65};
172  Int_t k = 0; // gmt
173  Double_t *rot = comb->GetRotationMatrix();
174  Double_t *tra = comb->GetTranslation();
175  const StThreeVectorD unit(1.,0.,0.);
176  const StThreeVectorD zero(0.,0.,0.);
177  StThreeVectorD normal(rot[2], rot[5], rot[8]);
178  StThreeVectorD middle(tra);
179  if (_debug) {
180  LOG_INFO << "middle:" << middle << "\tnormal:" << normal << endm;
181  }
182  comb->LocalToMaster(zero.xyz(),middle.xyz());
183  comb->LocalToMasterVect(unit.xyz(), normal.xyz());
184  if (_debug) {
185  LOG_INFO << "middle:" << middle << "\tnormal:" << normal << endm;
186  }
187  Double_t phiM = TMath::RadToDeg()*middle.phi();
188  if (_debug) {
189  LOG_INFO << "phiR = " << phiR << "\tphiM = " << phiM << endm;
190  }
191  Double_t dPhi = phiR - phiM;
192  if (dPhi > 360) dPhi -= 360;
193  if (dPhi < -360) dPhi += 360;
194  if (TMath::Abs(dPhi) > 15) continue;
195  if (_debug) {
196  LOG_INFO << "zR = " << xyzR.z() << "\tzM = " << tra[2] << endm;
197  }
198  if (TMath::Abs(xyzR.z() - tra[2]) > 20) continue;
199  Double_t sh = helixO.pathLength(middle, normal);
200  if (_debug) {
201  LOG_INFO << "StHelix sh " << sh
202  << "\t shR " << shR.first << "\t" << shR.second
203  << endm;
204  }
205  StThreeVectorD xyzG = helixO.at(sh); if (_debug) cout << "StHelix xyzG\t" << xyzG << endl;
206  StThreeVectorD dR = xyzR - xyzG; if (_debug) cout << "dR\t" << dR << " dist = " << dR.magnitude() << endl;
207  if (dR.magnitude() > 50) continue;
208  if (sh < -5e2 || sh > 5e2) continue;
209  if (_debug) {
210  StThreeVectorD dX = xyzG - helixO.at(0);
211  LOG_INFO << "Qi: " << Track->geometry()->charge()
212  << "\tQo: " << Track->outerGeometry()->charge()
213  << "\tdX " << dX << endm;
214  LOG_INFO << *dca << endm;
215  }
216  Double_t uvPred[3];
217  comb->MasterToLocal(xyzG.xyz(),uvPred);
218  TRVector xyzL(3,uvPred); if (_debug) cout << "StHelix xyzL\t" << xyzL << endl;
219  Double_t dirGPred[3] = {helixO.cx(sh),helixO.cy(sh),helixO.cz(sh)};
220  Double_t dxyzL[3];
221  comb->MasterToLocalVect(dirGPred,dxyzL);
222  Double_t tuvPred[2] = {dxyzL[1]/dxyzL[0], dxyzL[2]/dxyzL[0]};
223  if (_debug) {
224  LOG_INFO << "StHelix tU/tV = " << tuvPred[0] << "\t" << tuvPred[1] << endm;
225  }
226  if (TMath::Abs(uvPred[1]) > dx[k] + 1.0) continue;
227  if (TMath::Abs(uvPred[2]) > dz[k] + 1.0) continue;
228 
229  Bool_t mIsCrossingMembrain = kTRUE;
230  Bool_t mIsPrimary = kFALSE;
231  if(pTrack!=0) {
232  mIsPrimary = kTRUE;
233  }
234  StThreeVectorF firstPoint = Track->detectorInfo()->firstPoint();
235  StThreeVectorF lastPoint = Track->detectorInfo()->lastPoint();
236  if( (firstPoint.z()>0 && lastPoint.z()>0 && xyzG.z()>0) ||
237  (firstPoint.z()<0 && lastPoint.z()<0 && xyzG.z()<0) ) {
238  mIsCrossingMembrain = kFALSE;
239  }
240 
241  for (UInt_t l = 0; l < NoHits; l++) {
242  StHit *hit = hitvec[l];
243  if (hit) {
244  //if (hit->flag()>=4) continue;
245  //if (hit->flag()< 0) continu;
246  // LOG_INFO << "hitFlag=" << hit->flag() << endm;
247 
248  HitT *h = AddHitT();
249  h->SetHitLength(sh);
250  h->SetHitLength(stepR);
251  h->SetHitdR(dR.magnitude());
252  h->SetHitFlag(UInt_t(hit->flag()));
253  h->SetUVPred (uvPred[1],uvPred[2]);
254  h->SettUVPred(tuvPred[0],tuvPred[1]);
255  h->SetXyzG(xyzG.xyz());
256  h->SetDirG(dirGPred);
257  h->SetisPrimary(mIsPrimary);
258  h->SetisCrossingMembrain(mIsCrossingMembrain);
259  SetHitT(h, hit, comb, track);
260  NoHitPerTrack++;
261  h->SetHitPerTrack(NoHitPerTrack);
262  // SetHitT(h, hit, comb, track, &TPDeriv);
263  } // if (hit)
264  } // for (UInt_t l = 0; l < NoHits; l++)
265  } // while ((comb = (TGeoHMatrix *) next()))
266  nTotMatch += NoHitPerTrack;
267  } // for (UInt_t i=0; i<nnodes; i++)
268  fIsValid = kTRUE;
269  if (nTotMatch) iok = 0;
270  return iok;
271 }
272 
273 //_________________
274 TrackT *EventT::AddTrackT() {
275  // Add a new track to the list of tracks for this event.
276  // To avoid calling the very time consuming operator new for each track,
277  // the standard but not well know C++ operator "new with placement"
278  // is called. If tracks[i] is 0, a new TrackT object will be created
279  // otherwise the previous TrackT[i] will be overwritten.
280 
281  TClonesArray &tracks = *fTracks;
282  TrackT *track = new(tracks[fNtrack++]) TrackT();
283  //Save reference to last TrackT in the collection of Tracks
284  return track;
285 }
286 
287 //_________________
288 HitT *EventT::AddHitT() {
289  // Add a new hit to the list of hits for this event.
290  // To avoid calling the very time consuming operator new for each hit,
291  // the standard but not well know C++ operator "new with placement"
292  // is called. If hits[i] is 0, a new HitT object will be created
293  // otherwise the previous HitT[i] will be overwritten.
294 
295  TClonesArray &hits = *fHits;
296  HitT *hit = new(hits[fNhit++]) HitT();
297  //Save reference to last HitT in the collection of Hits
298  return hit;
299 }
300 
301 //_________________
302 void EventT::Clear(Option_t * /*option*/) {
303  fTracks->Clear("C"); //will also call TrackT::Clear
304  fHits->Clear("C"); //will also call HitT::Clear
305 }
306 
307 //_________________
308 void EventT::Reset(Option_t * /*option*/) {
309  // Static function to reset all static objects for this event
310  // fgTracks->Delete(option);
311 
312  delete fgTracks; fgTracks = 0;
313  delete fgHits; fgHits = 0;
314 }
315 
316 //_________________
317 void EventT::SetHeader(Int_t i, Int_t run, Int_t date, Double32_t field) {
318  fNtrack = 0;
319  fNhit = 0;
320  fEvtHdr.Set(i, run, date, field);
321 }
322 
323 //___________________
324 void EventT::Print(Option_t *opt) const {
325  LOG_INFO << "Run/EventT\t" << fEvtHdr.GetRun() << "/" << fEvtHdr.GetEvtNum() << "\tDate " << fEvtHdr.GetDate()
326  << "\tField " << fEvtHdr.GetField() << endm;
327  LOG_INFO << "Total no. tracks " << GetTotalNoTracks() << "\tRecorded tracks " << GetNtrack()
328  << "\tRecorded hits " << GetNhit() << endm;
329  TRVector vertex(3,GetVertex());
330  TRSymMatrix cov(3,GetCovMatrix());
331  LOG_INFO << "Primary vertex " << vertex << endm;
332  LOG_INFO << "Its cov. matrix " << cov << endm;
333  for (UInt_t i = 0; i < GetNtrack(); i++) {
334  LOG_INFO << i << "\t"; GetTrackT(i)->Print();
335  }
336  for (UInt_t i = 0; i < GetNhit(); i++) {
337  LOG_INFO << i << "\t"; GetHitT(i)->Print();
338  }
339 }
340 
341 //___________________
342 HitT *EventT::SetHitT(HitT *h, StHit *hit, TGeoHMatrix *comb, TrackT *track) {
343  UInt_t B = 0, L = 0, l = 0, W = 0, H = 0;
344  Int_t rdo = 0;
345  h->SetRDO(rdo);
346  if (hit->detector() == kGmtId) {
347  StGmtHit *ht = (StGmtHit *) hit;
348  B = ht->getModule();
349  h->SetId(B,L,l,W,H);
350  h->SetuvD(ht->getLocalY(), ht->getLocalX());
351  h->SetuvDError(ht->getErrorLocalY(), ht->getErrorLocalX());
352  h->SetSigma(ht->getSigmaY(), ht->getSigmaX());
353  h->SetSigmaError(ht->getErrorSigmaY(), ht->getErrorSigmaX());
354  h->SetAdc(ht->getAdcY(), ht->getAdcX());
355  h->SetAdcError(ht->getErrorAdcY(), ht->getErrorAdcX());
356  h->SetUsedInFit(hit->usedInFit());
357  } // if (hit->detector() == kGmtId)
358  StThreeVectorF position = hit->position();
359  Double_t xyzG[3] = {position.x(),position.y(),position.z()};
360  h->SetGC(xyzG[0],xyzG[1],xyzG[2]);
361  Double_t xyzL[3] = {0,0,0};
362  comb->MasterToLocal(xyzG,xyzL);
363  // if (TMath::Abs(xyzL[2]) > 0.1) continue;
364  Double_t uvw[3] = {h->GetU(),h->GetV(),0};
365  comb->LocalToMaster(uvw,xyzG);
366  h->Set(xyzG,uvw);
367  Double_t *rot = comb->GetRotationMatrix();
368  h->SetWG(rot[2],rot[5],rot[8]);
369  // Int_t IdH = GetIndexOfHitT(h);
370  Int_t IdH = fNhit - 1;
371  track->SetHitTId(IdH);
372  Double_t invpT = track->GetInvpT();
373  if (TMath::Abs(invpT) < 1e-7) invpT = 1e-7;
374  h->SetpT(1./invpT);
375  h->SetMom(track->GetMomentum());
376  h->SetWG(rot[2],rot[5],rot[8]);
377  TGeoHMatrix *rotL = (TGeoHMatrix *) RotMatrices()->FindObject(Form("WL%s",comb->GetName()+1));
378  Double_t xyzLadder[3] = {0,0,0};
379  if (rotL) {
380  rotL->LocalToMaster(uvw,xyzLadder);
381  h->SetL(xyzLadder[0],xyzLadder[1],xyzLadder[2]);
382  Double_t uvwP[3] = {h->GetPredU(),h->GetPredV(),0};
383  rotL->LocalToMaster(uvwP,xyzLadder);
384  h->SetXyzL(xyzLadder);
385  }
386  else {
387  LOG_INFO << Form("WL%s",comb->GetName()+1) << " has not been found" << endm;
388  h->SetL(xyzLadder[0],xyzLadder[1],xyzLadder[2]);
389  h->SetXyzL(xyzLadder);
390  }
391  return h;
392 }
393 
394 //___________________
395 void TrackT::Print(Option_t *opt) const {
396  LOG_INFO << "TrackT: InvpT " << fInvpT << "\tTanL " << fTanL
397  << "\tPhi " << fPhi << "\tRho " << fRho
398  << "\tNpoint " << fNpoint << "\tNsp " << fNsp << endm;
399  for (UInt_t i = 0; i < fNsp; i++) {
400  LOG_INFO << "\t" << fIdHitT[i];
401  }
402  LOG_INFO << endm;
403 }
404 
405 //___________________
406 void HitT::SetId(Int_t B, Int_t L, Int_t l, Int_t W, Int_t H) {
407  barrel = B; layer = L; ladder = l; wafer = W; hybrid = H;
408 }
409 
410 //___________________
411 void HitT::Print(Option_t *opt) const {
412  LOG_INFO << "HitT: Id " << Id << "\tpT = " << pT << "\tmomentum " << pMom << endm;
413  TRVector glob(3,&xG);
414  LOG_INFO << "Global :" << glob << endm;
415  LOG_INFO << "Local u/v/w " << u << "/ " << v << "/ " << w << endm;
416  LOG_INFO << "Prediction uP/vP " << uP << "/ " << vP << "\ttuP/tvP " << tuP << "/ " << tvP << endm;
417 }
418 
419 //___________________
420 void EventT::RestoreListOfRotations() {
421  if (fRotList) return;
422  if (! gDirectory) return;
423  fRotList = new THashList(100,0);
424  fRotList->SetOwner();
425  TIter nextkey(gDirectory->GetListOfKeys() );
426  TKey *key;
427  while ((key = (TKey*) nextkey())) {
428  TObject *obj = key->ReadObj();
429  if ( obj->IsA()->InheritsFrom( "TGeoHMatrix" ) ) {
430  fRotList->Add(obj);
431  }
432  }
433 }
434 
435 //___________________
436 void TBase::Loop(Int_t Nevents) {
437 
438 //#if 1
439 
440  struct PlotPar_t {
441  const Char_t *Name;
442  const Char_t *Title;
443  Int_t nx;
444  Int_t ny;
445  Double_t xmin;
446  Double_t xmax;
447  Double_t ymin;
448  Double_t ymax;
449  };
450  // plots for uP
451  const PlotPar_t plotUP = { "uP", "track u", 320, 3, -5., 5., 0., 3. };
452  // plots for u-uP
453  const PlotPar_t plotDu = { "Du", "Du before cut", 250, 3, -2., 2., 0., 3. };
454  // plots for du & dv
455  const PlotPar_t plotDuv = { "Du", "Du cuts", 200, 3, -1., 1., 0., 3. };
456 
457  TFile *fOut = new TFile(fOutFileName,"recreate");
458  TString Name;
459  TString Title;
460  TString uName;
461  TString uTitle;
462  enum {NM = 8}; // no. of modules
463  // B
464  TH1F *LocPlots[NM];
465  TH1F * uPlots[NM];
466 #ifdef Extra
467  TH1F *hpT = new TH1F( "Pt", "pt", 200, -2., 2.);
468  TH1F *hpM = new TH1F( "Ptot", "ptot", 200, 0., 5.);
469  TH1F *uAll = new TH1F("Uall", "ua", plotUP.nx, plotUP.xmin, plotUP.xmax);
470  TH1F * xCuts[NM];
471  TH1F * uCuts[NM];
472 #endif
473  TH1F * uPAll = new TH1F("UPall","uPall", plotUP.nx, plotUP.xmin, plotUP.xmax);
474  TH1F * duB[NM][2];
475  TH1F * dvB[NM];
476  TH1F * uCut = new TH1F("Ucut","uc", plotDu.nx, plotDu.xmin, plotDu.xmax);
477  TH1F * vCut = new TH1F("Vcut","vc", 200, -3., 3.);
478  TH2F * dMin = new TH2F("DMin","vumin",100,-0.75,0.75,100,-0.75,0.75);
479  TH1F * vMin = new TH1F("VMin","vmin", plotDuv.nx, plotDuv.xmin, plotDuv.xmax);
480  TH1F * uMin = new TH1F("UMin","umin", plotDuv.nx, plotDuv.xmin, plotDuv.xmax);
481  //TH1F * uMinC = new TH1F("UMinC","umC", plotDuv.nx, plotDuv.xmin, plotDuv.xmax);
482  memset(LocPlots,0,NM*sizeof(TH1F *));
483  memset( uPlots,0,NM*sizeof(TH1F *));
484  // Loop over gmt modules
485  for (int M = 0; M < NM; M++) {
486  uName = Form("UModule%i", M+1);
487  uTitle = Form("du for M%i", M+1);
488  duB[M][0] = new TH1F(uName, uTitle, plotDuv.nx, plotDuv.xmin, plotDuv.xmax );
489  uName = Form("VModule%i", M+1);
490  uTitle = Form("dv for M%i", M+1);
491  dvB[M] = new TH1F(uName, uTitle, plotDuv.nx, plotDuv.xmin, plotDuv.xmax );
492  uName = Form("UModule%iVcut", M+1);
493  uTitle = Form("du for M%i after Vcut", M+1);
494  duB[M][1] = new TH1F(uName, uTitle, plotDuv.nx, plotDuv.xmin, plotDuv.xmax );
495  Int_t module = M;
496  uName = plotUP.Name;
497  uName += Form("M%i", module);
498  uTitle = Form("uP for Module %i", module);
499  uPlots[module] = new TH1F(uName, uTitle, plotUP.nx, plotUP.xmin, plotUP.xmax );
500  uName = Form("%sM%i", plotDu.Name, module);
501  uTitle = Form("u-uP for M %i", module);
502  uName = Form("%sxM%i", plotDu.Name, module);
503  uTitle = Form("u-uP corr M %i", module);
504 #ifdef Extra
505  xCuts[module] = new TH1F(uName, uTitle, plotDu.nx, plotDu.xmin, plotDu.xmax );
506  uCuts[module] = new TH1F(uName, uTitle, plotDu.nx, plotDu.xmin, plotDu.xmax );
507 #endif
508  } // for (int M = 0; M < NM; M++)
509 
510  Long64_t nentries = fChain->GetEntriesFast();
511  if (Nevents > 0 && nentries > Nevents) {
512  nentries = Nevents;
513  }
514 
515  Long64_t nbytes = 0, nb = 0;
516  Int_t TreeNo = -1;
517  TString currentFile("");
518 
519  // Loop over events in tree
520  for (Long64_t jentry=0; jentry<nentries;jentry++) {
521 
522  Long64_t ientry = LoadTree(jentry);
523  if (ientry < 0) break;
524  nb = fChain->GetEntry(jentry); nbytes += nb;
525  if (! jentry%1000 || TreeNo != fChain->GetTreeNumber()) {
526  LOG_INFO << "Read event \t" << jentry
527  << " so far, switch to file " << fChain->GetCurrentFile()->GetName()
528  << endm;
529  LOG_INFO << " current TreeNo: " << TreeNo
530  << " new TreeNo: " << fChain->GetTreeNumber() << endm;
531  TreeNo = fChain->GetTreeNumber();
532  } // for (Long64_t jentry=0; jentry<nentries;jentry++)
533 
534  // if (VertexZCut > 0 && TMath::Abs(fVertex[2]) > VertexZCut) continue;
535  UInt_t Ntrack = fEvent->GetTracks()->GetEntriesFast();
536  // int k_used[100000] = {0};
537 
538  // Loop over tracks
539  for (UInt_t trk = 0; trk < Ntrack; trk++) {
540  TrackT *track = (TrackT *)fEvent->GetTracks()->UncheckedAt(trk);
541  if (! track) continue;
542 #ifdef Extra
543  Int_t Npoints = track->GetNpoint();
544  if (minNoFitPoints > 0 && Npoints%100 < minNoFitPoints) continue;
545  if (UseSsd && Npoints < 1000) continue;
546  if (UseSvt && Npoints < 100) continue;
547 #endif
548  Int_t Nsp = track->GetN();
549  double dvmin = 1000.;
550  double dumin = 1000.;
551  //int kmin{0};
552  // Loop over hits
553  for (Int_t hit = 0; hit < Nsp; hit++) {
554  Int_t k = track->GetHitTId(hit) - 1;
555  // assert(k>=0);
556  HitT *hitT = (HitT *) fEvent->GetHitT(k);
557  if ( k < 0 ) LOG_INFO <<" k <0:"<<k<<" hit="<<hit<<" Nsp="<<Nsp<< endm;
558  if ( k < 0 ) continue;
559  Int_t module = hitT->Barrel();
560  Double32_t u = hitT->GetU();
561  Double32_t v = hitT->GetV();
562  Double32_t uP = hitT->GetPredU();
563  Double32_t vP = hitT->GetPredV();
564  Double32_t du = u - uP;
565  Double32_t dv = v - vP;
566 #ifdef Extra
567  hpT->Fill(hitT->GetXyzPGl());
568  hpM->Fill(hitT->GetpMom());
569  if (TMath::Abs(hitT->GetpT()) < 0.2) continue;
570 #endif
571 
572  if ( TMath::Abs(dv) < TMath::Abs(dvmin) ) {
573  dvmin = dv;
574  dumin= du;
575  //kmin = k;
576  }
577  uPAll->Fill( uP );
578  uPlots[module]->Fill( uP );
579 
580  duB[module][0]->Fill(du);
581  dvB[module]->Fill(dv);
582  vCut->Fill(dv);
583  // if (TMath::Abs(dv) > rCut ) continue;
584  uCut->Fill(du);
585 
586  duB[module][1]->Fill(du);
587  // if (TMath::Abs(du) > 2.*rCut) continue;
588  } //hits loop
589 
590  if (TMath::Abs(dvmin) < 1000.) {
591  dMin->Fill(dvmin,dumin);
592  vMin->Fill(dvmin);
593  uMin->Fill(dumin);
594  // if (TMath::Abs(dvmin) < rCut ) uMinC->Fill(dumin);
595  } // if (TMath::Abs(dvmin) < 1000.)
596 
597  } //track loop
598  } //jentry loop (event loop)
599  fOut->Write();
600 
601 //#endif // #if 1
602 
603 } //end of Loop()
Holds collections of GMT data.
Float_t getErrorLocalX() const
Local X coordinate error.
Definition: StGmtHit.h:60
Float_t getErrorSigmaX() const
Position error in X.
Definition: StGmtHit.h:68
Float_t getErrorAdcY() const
ADC error in Y.
Definition: StGmtHit.h:56
StGmtHitCollection * getHitCollection(unsigned short moduleIdx)
Pointer to the GMT hit collection in the i-th module.
Float_t getSigmaY() const
Position in Y.
Definition: StGmtHit.h:70
Definition: StHit.h:125
Float_t getAdcX() const
ADC in X.
Definition: StGmtHit.h:50
Definition: TrackT.h:15
Float_t getLocalX() const
Local X coordinate.
Definition: StGmtHit.h:58
Definition: HitT.h:16
EventT * fEvent
current Tree number in a TChain
Definition: EventT.h:106
Float_t getLocalY() const
Local Y coordinate.
Definition: StGmtHit.h:62
Int_t getModule() const
Module.
Definition: StGmtHit.h:48
Holds collection of GMT hits in the module.
double cx(double s) const
pointing vector of helix at point s
Definition: StHelix.hh:203
Holds data for the hit in GMT.
Definition: StGmtHit.h:23
pair< double, double > pathLength(double r) const
path length at given r (cylindrical r)
Definition: StHelix.cc:351
Float_t getSigmaX() const
Position in X.
Definition: StGmtHit.h:66
Float_t getErrorLocalY() const
Local Y coordinate error.
Definition: StGmtHit.h:64
Float_t getAdcY() const
ADC in Y.
Definition: StGmtHit.h:54
StSPtrVecGmtHit & getHitVec()
Vector of hits that belong to the module.
C++ STL includes.
Definition: AgUStep.h:47
Float_t getErrorAdcX() const
ADC error in X.
Definition: StGmtHit.h:52
Float_t getErrorSigmaY() const
Position error in X.
Definition: StGmtHit.h:72