StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StiForwardTrackMaker.cxx
1 // -- Author : Victor Perevoztchikov
2 //
3 // $Id: StiForwardTrackMaker.cxx,v 1.10 2007/07/12 19:27:21 fisyak Exp $
4 #include "TFile.h"
5 #include "StMessMgr.h"
6 
7 #include "StiForwardTrackMaker.h"
8 #include "StEvent.h"
9 #include "StPrimaryVertex.h"
10 
11 #include "Sti/StiToolkit.h"
12 #include "Sti/StiKalmanTrack.h"
13 #include "Sti/StiKalmanTrackNode.h"
14 #include "StiMaker/StiMaker.h"
15 #include "Sti/StiTrackContainer.h"
16 
17 #include "Sti/StiHitContainer.h"
18 #include "Sti/StiLocalTrackSeedFinder.h"
19 #include "Sti/StiHit.h"
20 #include "Sti/StiSortedHitIterator.h"
21 #include "Sti/StiDetectorContainer.h"
22 #include "Sti/StiKalmanTrackFinder.h"
23 
24 #include "StDetectorId.h"
25 #include "StHit.h"
26 #include "TFile.h"
27 #define xG(t) (t->x_g())
28 #define yG(t) (t->y_g())
29 #define zG(t) (t->z_g())
30 #define rxyG(t) sqrt(xG(t)*xG(t) + yG(t)*yG(t))
31 #define rG(t) sqrt(xG(t)*xG(t) + yG(t)*yG(t) + zG(t)*zG(t))
32 #define xL(t) (t->getX())
33 #define yL(t) (t->getY())
34 #define zL(t) (t->getZ())
35 #define ezL(t) sqrt(t->getCzz())
36 #define eyL(t) sqrt(t->getCyy())
37 
38 
39 ClassImp(StiForwardTrackMaker)
40 
41 //_____________________________________________________________________________
42 StiForwardTrackMaker::StiForwardTrackMaker(const char *name):StMaker(name){
43  mToolkit=0;
44  mTotEve=0;
45  mTotMatchedSeeds=0;
46  mTotMatchedTracks=0;
47  HList=0;
48  memset(hA,0,sizeof(hA));
49  mMaxTrkDcaRxy = 3.0; // cm
50  mMaxZdca = 4; // cm
51 
52  mMinEta = 0.9;
53  mAllHits = 0;
54  mForwardHits = 0;
55  mTrackSeeds = 0;
56  mSeedGenerator = 0;
57  mTrackFinder = 0;
58 }
59 
60 
61 //_____________________________________________________________________________
62 StiForwardTrackMaker::~StiForwardTrackMaker(){
63  //
64 }
65 
66 
67 //_____________________________________________________________________________
70 
71  //get a pointer to StiMaker:
72  StiMaker* sti = (StiMaker*)StMaker::GetChain()->GetMaker("Sti");
73  if(sti==0) {
74  gMessMgr->Warning() <<GetName()<<"no STi Maker, it is fatal"<<endm;
75  return kStErr;
76  }
77 
78  //get pointer to Sti toolkit
79  mToolkit = sti->getToolkit();
80  assert(mToolkit); //internal error of Sti
81 
82  //get pointer to full hit container
83  mAllHits = mToolkit->getHitContainer();
84  assert(mAllHits);
85 
86  //instantiate the forward hit container
87  mForwardHits = new StiHitContainer("mForwardHits","forward hits",mToolkit->getHitFactory());
88 
89  //instantiate the track seeds container
90  mTrackSeeds = new StiTrackContainer("mTrackSeeds","forward track seeds");
91 
92  //instantiate the seed generator
93  mSeedGenerator = new StiLocalTrackSeedFinder("mSeedGenerator","builds track seeds",mToolkit->getTrackFactory(), mForwardHits,mToolkit->getDetectorContainer());
94 
95  //initialize the Kalman track finder
96  mTrackFinder = new StiKalmanTrackFinder(mToolkit);
97  //mTrackFinder->initialize();
98 
99  HList=new TObjArray(0);
100  initHisto();
101  HList->ls();
102 
103  gMessMgr->Info() << GetName()
104  <<"::Cuts"
105  <<"\n MaxTrkDcaRxy/cm = "<< mMaxTrkDcaRxy
106  <<"\n mMaxZdca = "<< mMaxZdca
107  <<"\n mMinEta = "<< mMinEta
108  <<endm;
109 
110  return StMaker::Init();
111 }
112 
113 
114 //_____________________________________________________________________________
117  return MakeAfterSti();
118 }
119 
120 //_____________________________________________________________________________
123 {
124  if(mToolkit==0) {
125  gMessMgr->Warning()<<GetName()<<"no Sti tool kit "<<GetName()<<" is OFF"<<endm;
126  return kStErr;
127  }
128  else if(mAllHits==0) {
129  gMessMgr->Warning() <<"no Sti hits present: "<<GetName()<<" is OFF"<<endm;
130  return kStErr;
131  }
132 
133 
134  int nV=vertL.size();
135  StEvent *event = (StEvent *) GetInputDS("StEvent");
136 
137  gMessMgr->Info() << "\n JJJ1 "<<GetName()<<":: MakeInSti(), eveID="<<event->id()<<" nVert="<<nV<<endm;
138 
139  //print the position of each of the vertices, fill histo with zVertex
140  int iv;
141  for(iv=0;iv<nV;iv++) {
142  VertexV V=vertL[iv];
143  printf("iv=%d Vz=%.2f +/-%.2f \n",iv,V.z,V.ez );
144  hA[1]->Fill(V.z);
145  }
146 
147  getForwardHits(mAllHits, mForwardHits, mToolkit->getDetectorContainer(), mMinEta);
148  cout<<"We started with " << mAllHits->size() << " hits and ended with " << mForwardHits->size() << " in the forward region" << endl;
149 
150  buildTrackSeeds(mForwardHits, mTrackSeeds, mSeedGenerator);
151 
152  //use Kalman to extend tracks, build covariance matrix, and add more hits
153  extendTracks();
154 
155 
156  //make some histograms
157  StiSortedHitIterator it = StiSortedHitIterator(mForwardHits, mToolkit->getDetectorContainer()->begin(), mToolkit->getDetectorContainer()->end());
158  StiHit* hit = 0;
159  float phi;
160  while(it!=StiSortedHitIterator())
161  {
162  hit = &(*it);
163 
164  //eta
165  hA[2]->Fill( -log( rxyG(hit)/(zG(hit)+rG(hit)) ) ); //eta of hit
166 
167  //phi (is this calculation correct?)
168  phi = atan(yG(hit)/xG(hit));
169  if(xG(hit)<0 && yG(hit)<0) phi = phi - 3.14;
170  else if(xG(hit)<0) phi = phi + 3.14;
171  hA[3]->Fill(phi);
172 
173  it++;
174  }
175 
176 
177  //get the Sti track container...
178  StiTrackContainer* tracks = mToolkit->getTrackContainer();
179  if(tracks==0) {
180  gMessMgr->Warning() <<"no STi tracks , skip eve"<<endm;
181  return kStErr ;
182  }
183 
184  cout<<"There are "<<tracks->getTrackCount(0)<<" original tracks and "<<mTrackSeeds->getTrackCount(0)<<" forward track seeds"<<endl;
185 
186  matchVertex(tracks,vertL,mMaxZdca,nV,hA[9],mTotMatchedTracks);
187 
188  matchVertex(mTrackSeeds,vertL,mMaxZdca,nV,hA[8],mTotMatchedSeeds);
189 
190  return kStOK;
191 }
192 
193 //_____________________________________________________________________________
196 
197  StEvent *event = (StEvent *) GetInputDS("StEvent");
198  assert(event);
199  eveID=event->id();
200  mTotEve++;
201  int nV=event->numberOfPrimaryVertices();
202 
203  gMessMgr->Info() << "\n JJJ2 "<<GetName()<<"MakeAfterSti(), START nEve="<<mTotEve<<" eveID="<<eveID<< " nPrimVertex="<<nV<< endm;
204 
205  hA[0]->Fill(nV);
206 
207  if(nV<=0) {
208  gMessMgr->Info() << GetName()<<" event dropped, no vertex found"<<endm;
209  return kStOK;
210  }
211 
212 
213  return kStOK;
214 }
215 
216 //_____________________________________________________________________________
217 // called before each event
218 void
219 StiForwardTrackMaker::Clear(const char* opt){
220  gMessMgr->Info() <<GetName()<<"::Clear()"<< endm;
221  eveID=0;
222  vertL.clear();
223  mForwardHits->clear();
224  mTrackSeeds->clear();
225  mSeedGenerator->reset();
226 }
227 
228 //=======================================================
229 //=======================================================
230 void
231 StiForwardTrackMaker::initHisto() {
232  assert(HList);
233  hA[0]=new TH1F("nV","No. of vertices per eve",20,-0.5,19.5);
234  hA[1]=new TH1F("zV","reconstructed vertices (any); Z (cm)",200,-200,200);
235 
236  hA[2] = new TH1F("etaNotUsed","eta of unused forward hits",20,0.8,2.1);
237  hA[3] = new TH1F("phiNotUsed","phi of unused forward hits",20,-3.14,3.14);
238  hA[4] = new TH1F("etaUsed","eta of hits used by previous tracker",20,0.8,2.1);
239  hA[5] = new TH1F("phiUsed","phi of hits used by previous tracker",20,-3.14,3.14);
240  hA[6] = new TH1F("seedHits","number of seed hits used for each track",20,1,20);
241  hA[7] = new TH1F("trackHits1","total # of hits/track after inward Kalman",20,1,20);
242  hA[10] = new TH1F("trackHits2","total # of hits/track after outward Kalman",20,1,20);
243  hA[8] = new TH1F("deltaZkalman","zDCA-zVertex for tracks after Kalman",200,-200,200);
244  hA[9] = new TH1F("deltaZtracks","zDCA-zVertex for old tracks",200,-200,200);
245 
246  int i;
247  for(i=0;i<mxHA; i++) if(hA[i]) HList->Add(hA[i]);
248 }
249 
250 //-------------------------------------------------
251 //-------------------------------------------------
252 void
253 StiForwardTrackMaker::saveHisto(TString fname){
254  TString outName=fname+".hist.root";
255  TFile f( outName,"recreate");
256  assert(f.IsOpen());
257  printf("%d histos are written to '%s' ...\n",HList->GetEntries(),outName.Data());
258  HList->ls();
259  HList->Write();
260  f.Close();
261 }
262 
263 //-------------------------------------------------
264 //-------------------------------------------------
265 Int_t
267  saveHisto(GetName());
268  cout<<"Total matched track seeds: "<<mTotMatchedSeeds<<" and total matched tracks: "<<mTotMatchedTracks<<endl;
269  return kStOK;
270 }
271 
272 //-------------------------------------------------
273 //-------------------------------------------------
274 
275 void StiForwardTrackMaker::getForwardHits(StiHitContainer* allHits, StiHitContainer* forwardHits, StiDetectorContainer* detector, double minEta)
276 {
277  StiSortedHitIterator it = StiSortedHitIterator(allHits, detector->begin(), detector->end());
278  StiHit* hit = 0;
279 
280  while(it!=StiSortedHitIterator())
281  {
282  hit = &(*it);
283  StDetectorId id = (static_cast <const StHit*>(hit->stHit()))->detector(); //restricting to TPC hits; doesn't seem to do anything
284  if( rxyG(hit)/(zG(hit)+rG(hit)) < exp(-minEta) && id==kTpcId)
285  {
286  if(hit->timesUsed()==0)forwardHits->add(hit);
287  else //this is a previously used hit, just put it in a histogram
288  {
289  hA[4]->Fill( -log( rxyG(hit)/(zG(hit)+rG(hit)) ) ); //eta of hit
290 
291  //phi (is this calculation correct?)
292  float phi = atan(yG(hit)/xG(hit));
293  if(xG(hit)<0 && yG(hit)<0) phi = phi - 3.14;
294  else if(xG(hit)<0) phi = phi + 3.14;
295  hA[5]->Fill(phi);
296  }
297  }
298  it++;
299  }
300 
301  forwardHits->sortHits();
302 }
303 
304 
305 void StiForwardTrackMaker::buildTrackSeeds(const StiHitContainer* forwardHits, StiTrackContainer* trackSeeds, StiLocalTrackSeedFinder* seedGenerator)
306 {
307  StiTrack* track = seedGenerator->findTrack();
308 
309  while(track != 0)
310  {
311  hA[6]->Fill(track->getPointCount());
312  trackSeeds->add(track);
313  track = seedGenerator->findTrack();
314  }
315 }
316 
317 void StiForwardTrackMaker::extendTracks()
318 {
319  mTrackFinder->initialize();
320  //mTrackFinder->clear();
321 
322  float zDCA, ezDCA, RxyDCA;
323  //int goodTracks = 0;
324 
325  for(StiTrackContainer::const_iterator it=mTrackSeeds->begin(); it!=mTrackSeeds->end(); it++)
326  {
327  while(mTrackFinder->find(*it,0,2)){}; //extend inward
328  hA[7]->Fill((*it)->getPointCount());
329  while(mTrackFinder->find(*it,1,0)){}; //extend outward
330  hA[10]->Fill((*it)->getPointCount());
331 
332  if(examineTrackDca(static_cast<StiKalmanTrack*>(*it), zDCA, ezDCA, RxyDCA)){
333  // goodTracks++;
334  cout << "This track passed the cut "<<**it<<endl;
335  }
336  }
337 
338  //cout<<goodTracks<<" of the forward track seeds passed the examineTrackDca() cut"<<endl;
339 }
340 
341 
342 
343 
344 void StiForwardTrackMaker::matchVertex(StiTrackContainer* tracks, vector <VertexV> &vertL, double &mMaxZdca, int &nV, TH1* h, int &totalMatched)
345 {
346  int nAll=0, nAny=0, nTry=0, nAcc=0;
347  for (StiTrackContainer::const_iterator it=(*tracks).begin(); it!=(*tracks).end(); ++it) {
348  const StiKalmanTrack* track = static_cast<StiKalmanTrack*>(*it);
349  nAll++;
350  if(track->getFlag()!=true) nAny++; // don't drop bad tracks for now, all track seeds currently flagged as bad
351  //cout<<"\n#a kalTrack: nTr="<<nAny<<" flag="<<track->getFlag()<<" nFitP="<<track->getFitPointCount()<<" is Prim="<<track->isPrimary()<<" pT="<<track->getPt()<<endl;
352  // cout<<"#b kalTrack: pT="<<track->getPt()<<endl;
353  //cout<<"#b kalTrack:"<<*track<<endl;
354 
355  float zDca, ezDca, rxyDca;
356  if(!examineTrackDca(track, zDca, ezDca, rxyDca)) continue;
357  nTry++;
358  //cout<<"#c nTry="<<nTry<<" zDca="<<zDca<<endl;
359 
360  //.... match track to a vertex
361  for(int iv=0;iv<nV;iv++) {
362  VertexV V=vertL[iv];
363  h->Fill(zDca-V.z);
364  if( fabs(zDca-V.z) > mMaxZdca )continue;
365  nAcc++;
366  // cout<<" tr Matched , dZ="<<zDca-V.z<<" nAcc="<<nAcc<<endl;
367  break;
368  }
369  }
370 
371  gMessMgr->Info() << "\n"<<GetName()<<" We have "<<nAll<<" tracks, found "<<nAny<<" flagged tracks, try ZDca for "<<nTry<<", match to vertex:"<<nAcc<<"\n"<<endm;
372  totalMatched = nAcc + totalMatched;
373 }
374 
375 bool
376 StiForwardTrackMaker::examineTrackDca(const StiKalmanTrack*track,
377  float &zDca, float &ezDca, float &rxyDca){
378 
379  // cout <<"#e track->getPseudoRapidity()="<<track->getPseudoRapidity()<<" track->getFitPointCount()="<<track->getFitPointCount()<<endl;
380 
381  // .......... test DCA to beam .............
382  StiKalmanTrack track1=*track; // clone track
383  StiKalmanTrackNode* bmNode=track1.extrapolateToBeam();
384  if(bmNode==0) {
385  // cout<<"#a @beam DCA NULL"<<endl;
386  return false ;
387  }
388 
389  float rxy=rxyG(bmNode);
390  //cout<<"#e @beam global DCA x:"<< bmNode->x_g()<<" y:"<< bmNode->y_g()<<" z:"<< bmNode->z_g()<<" Rxy="<< rxy <<endl;
391  if(rxy>mMaxTrkDcaRxy) return false;
392 
393  //cout<<"#e inBeam |P|="<<bmNode->getP()<<" pT="<<bmNode->getPt()<<" local x="<<xL(bmNode)<<" y="<<yL(bmNode)<<" +/- "<<eyL(bmNode)<<" z="<<zL(bmNode)<<" +/- "<<ezL(bmNode)<<endl;
394 
395  zDca=zG(bmNode);
396  ezDca=ezL(bmNode);
397  rxyDca=rxy;
398 
399  return true;
400 }
401 
402 //-------------------------------------------------
403 //-------------------------------------------------
404 void
405 StiForwardTrackMaker::addVertex(float z, float ez){
406  VertexV V;
407  V.z=z;
408  V.ez=ez;
409  vertL.push_back(V);
410 }
411 
412 // $Log: StiForwardTrackMaker.cxx,v $
413 // Revision 1.10 2007/07/12 19:27:21 fisyak
414 // Add includes for TMath for ROOT 5.16
415 //
416 // Revision 1.9 2006/10/17 18:47:49 fisyak
417 // Remove reference to dead classes
418 //
419 // Revision 1.8 2006/08/06 23:29:44 perev
420 // ROOT5 corrs
421 //
422 // Revision 1.7 2005/09/26 14:03:54 kocolosk
423 // added Kalman propagation for forward track seeds
424 //
425 // Revision 1.6 2005/09/21 15:37:06 kocolosk
426 // modified matchVertex() so that |deltaZ| could be plotted separately for old, new tracks
427 //
428 // Revision 1.5 2005/09/14 14:15:06 kocolosk
429 // restrict to TPC hits, added histos, moved code to matchVertex() function
430 //
431 // Revision 1.4 2005/09/13 14:24:15 kocolosk
432 // implemented getForwardHits() and buildTrackSeeds() functions
433 //
434 // Revision 1.3 2005/09/12 21:08:21 balewski
435 // split Make to InSti and AfterSti
436 //
437 // Revision 1.2 2005/09/09 15:55:00 balewski
438 // prototype with hardcoded hacks
439 //
440 // Revision 1.1 2005/09/08 21:42:03 balewski
441 // star
442 //
443 
444 
445 
446 
447 
448 
449 
450 
virtual void initialize()
Initialize the finder.
virtual void add(StiHit *)
Definition of Kalman Track.
bool find(StiTrack *track, int direction, double rmin=0)
Find/extend the given track, in the given direction.
Definition: StHit.h:125
Abstract definition of a Track.
Definition: StiTrack.h:59
virtual Int_t Make()
Make - this method is called in loop for each event.
Definition: StiHit.h:51
Int_t MakeInSti()
Inside Sti.
virtual void reset()
Reset the tracker.
virtual void clear()
virtual void sortHits()
Definition of Kalman Track.
virtual Int_t Init()
Init - is a first method the top level StChain calls to initialize all its makers.
StiKalmanTrackNode * extrapolateToBeam()
StiTrack * findTrack(double rMin=0)
virtual unsigned int size() const
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:40
Abstract interface for a STI toolkit.
int getTrackCount(Filter< StiTrack > *filter) const
void add(StiTrack *track)
Add the given track to the container.
Definition: Stypes.h:44
Int_t MakeAfterSti()
Inside Sti.