StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StMiniMcMaker.cxx
Go to the documentation of this file.
1 
12 #include "StMiniMcMaker.h"
13 #include "TFile.h"
14 #include "TTree.h"
15 #include "TBranch.h"
16 #include "TMath.h"
17 #include "TRandom.h"
18 
19 
20 
21 #include "Stiostream.h"
22 #include <assert.h>
23 
24 #include "StMessMgr.h"
25 #include "PhysicalConstants.h"
26 #include "StPhysicalHelix.hh"
27 #include "SystemOfUnits.h"
28 #include "StIOMaker/StIOMaker.h"
29 #include "StParticleDefinition.hh"
30 #include "StMatrixF.hh"
31 #include "StChainOpt.h"
32 /*
33 #include "StPionPlus.hh"
34 #include "StPionMinus.hh"
35 #include "StProton.hh"
36 #include "StKaonMinus.hh"
37 #include "StKaonPlus.hh"
38 #include "StAntiProton.hh"
39 #include "StDeuteron.hh"
40 #include "StElectron.hh"
41 #include "StPositron.hh"
42 */
43 
44 #include "StEventTypes.h"
45 #include "StMcEventTypes.hh"
46 #include "StMcEvent.hh"
47 #include "StuRefMult.hh"
48 #include "StTpcDedxPidAlgorithm.h"
49 #include "StuProbabilityPidAlgorithm.h"
50 
51 #include "StEmcUtil/geometry/StEmcGeom.h"
52 #include "StEmcUtil/projection/StEmcPosition.h"
53 
54 #include "StMiniMcHelper.h"
55 
56 #include <vector>
57 #include <map>
58 #include "StTrack.h"
59 #include "StGlobalTrack.h"
60 #include "StDcaGeometry.h"
61 #include "StContainers.h"
62 #include "StTrackDetectorInfo.h"
63 #include "StEnumerations.h"
64 #include "StHit.h"
65 #include "THelixTrack.h"
66 
67 static int StMiniMcMakerErrorCount=0;
68 
69 //______________________________________________________________________________
70 //helper funtion prototypes
71 void dominatrackInfo(const StTrack*, int& dominatrackKey, short&, float&);
72 
73 //______________________________________________________________________________
74 // increasing order
75 bool hitCmp(const StTpcHit* p1,const StTpcHit* p2)
76 {
77  return (p1->position().perp()<p2->position().perp());
78 }
79 //______________________________________________________________________________
80 // decreasing energy
81 bool StMcCalorimeterHitCompdE(const StMcCalorimeterHit* const &lhs,const StMcCalorimeterHit* const &rhs) {
82  return (lhs->dE() > rhs->dE());
83 }
84 //______________________________________________________________________________
85 // decreasing energy
86 bool StEmcRawHitCompEne(const StEmcRawHit* lhs, const StEmcRawHit* rhs) {
87  return (lhs->energy() > rhs->energy());
88 }
89 
90 //______________________________________________________________________________
91 float scaleFactor(double Eta, int hitType=0) {
92 // hitType = 0 - towers
93 // hitType = 1 - pre shower
94 // hitType = 2 - shower max Eta
95 // hitType = 3 - shower max Phi
96  float P0[]={14.69,559.7,0.1185e6,0.1260e6};
97  float P1[]={-0.1022,-109.9,-0.3292e5,-0.1395e5};
98  float P2[]={0.7484,-97.81,0.3113e5,0.1971e5};
99 
100  float x=fabs(Eta);
101  return P0[hitType]+P1[hitType]*x+P2[hitType]*x*x;
102 }
103 
104 ClassImp(StMiniMcMaker)
105 
106 //---------CONSTRUCTORS, ETC--------
107 
108 //______________________________________________________________________________
109 StMiniMcMaker::StMiniMcMaker(const Char_t *name, const Char_t *title)
110  :
111  StMaker(name,title),
112  mMiniMcEvent(0),
113  mIOMaker(0),
114  mMiniMcTree(0),
115  mMiniMcDST(0),
116  mInFileName(),
117  mInFilePrefix(),
118  mOutFileName(),
119  mOutDir("./"),
120  mParameterFileName(),
121  mRcEvent(0),
122  mMcEvent(0),
123  mRcHitMap(0),
124  mRcTrackMap(0),
125  mMcTrackMap(0),
126  mRcVertexPos(0),
127  mMcVertexPos(0),
128  mTpcDedxAlgo(0),
129  mPidAlgo(0),
130  mEmcIndex(4801),
131  mGhost(kTRUE),
132  mMinPt(0),mMaxPt(99999),
133  mNSplit(0),mNRc(0),mNGhost(0),mNContam(0),
134  mNMatched(0),mNMatGlob(0), mMainVtx(-1)
135 
136 {
137 
138  // mParameterFileName = "/auto/data05/snelling/analysis/cvs/PIDTable.root";
139  // mParameterFileName = "/auto/pdsfdv05/starhipt/cbum/PIDTable.root";
140 
141 }
142 
143 //______________________________________________________________________________
144 StMiniMcMaker::~StMiniMcMaker()
145 {
146  /* */
147 }
148 
149 //
150 //---------- StMaker methods ----------
151 //
152 
153 /*
154  standard Clear()
155  */
156 //______________________________________________________________________________
157 void StMiniMcMaker::Clear(Option_t* opt)
158 {
159  StMaker::Clear();
160 }
161 
162 /*
163  called at the end
164 */
165 //______________________________________________________________________________
167 {
168  cout << "###StMiniMcMaker::Finish()" << endl;
169 
170  closeFile();
171 
172  cout << "\treconstr. = " << mNRc << endl
173  << "\tmatched = " << mNMatched << endl
174  << "\tsplit = " << mNSplit << endl
175  << "\tcontam. = " << mNContam << endl
176  << "\tghost = " << mNGhost << endl
177  << "\tmat global = " << mNMatGlob << endl;
178 
179  if (Debug()) cout << "deleting mMiniMcEvent" << endl;
180  delete mMiniMcEvent; mMiniMcEvent = 0;
181 // for some reason, the tree doesn't like to get deleted here...
182 // if (Debug()) cout << "deleting mMiniMcTree" << endl;
183 // delete mMiniMcTree;
184 // mMiniMcTree = 0;
185  if (Debug()) cout << "deleting mMiniMcDST" << endl;
186  delete mMiniMcDST;
187  mMiniMcDST = 0;
188 
189  return StMaker::Finish();
190 }
191 /*
192 
193  */
194 
195 //______________________________________________________________________________
196 Int_t StMiniMcMaker::InitRun(int runID)
197 {
198  cout << "###StMiniMcMaker::InitRun()" << endl;
199 
200  cout << "\tpt cut : " << mMinPt << " , " << mMaxPt << endl;
201  Int_t stat=0;
202 
203  //
204  // instantiate the event object here (embedding or simulation?)
205  //
206  if (!mMiniMcEvent) {
207  if(Debug()) cout << "\tStMiniMcMaker::InitRun Creating StMiniMcEvent..." << endl;
208  mMiniMcEvent = new StMiniMcEvent();
209  }
210  else {
211  if (Debug()) cout << "\tStMiniMcMaker::InitRun StMiniMcEvent Already created" << endl;
212  }
213  if(mGhost) {
214  cout << "\tGhost loop is turned on." << endl;
215  }
216 
217  //
218  // init the tpc dedx algo once
219  //
220  if (!mTpcDedxAlgo)
221  mTpcDedxAlgo = new StTpcDedxPidAlgorithm;
222 
223  //
224  // create file, trees, etc.
225  //
226  stat = openFile();
227 
228  return stat + StMaker::InitRun(runID);
229 
230 }
231 
232 //______________________________________________________________________________
233 Int_t StMiniMcMaker::Init()
234 {
235  //Moved everything important to InitRun(int)
236  cout << "###StMiniMcMaker::Init()" << endl;
237 
238  cout << "\tpt cut : " << mMinPt << " , " << mMaxPt << endl;
239 
240  if (mInFileName == "") {
241  const StChainOpt *chain = GetChainOpt();
242  assert(chain);
243  mInFileName = chain->GetFileOut();
244  }
245  return StMaker::Init();
246 }
247 
248 /*
249  Make called every event
250  */
251 
252 //______________________________________________________________________________
254 {
255  if(Debug()) cout << "###StMiniMcMaker::Make()" << endl;
256 
257 
258  //
259  // initialize StEvent, StMcEvent, and StAssociationMaker
260  //
261  mRcEvent = (StEvent*) GetDataSet("StEvent");
262  if(!mRcEvent) return kStOk; // last event apparently
263  mMcEvent = (StMcEvent*) GetDataSet("StMcEvent");
264  if(!mMcEvent) return kStErr;
265 
266  //
267  // association
268  //
269  Bool_t assOk = initAssociation();
270  if(!assOk) {
271  gMessMgr->Warning() << "Association problems " << endm;
272 // return kStErr;
273  }
274  //
275  // vertex
276  //
277  Bool_t vtxOk = initVertex();
278  if(!vtxOk) {
279  cout << "\t\t----No primary vertex---- " << endl;
280  return kStOk;
281  }
282 
283  // fill emc Helper array to index StEmcRawHits by
284  // SoftId.
285  buildEmcIndexArray();
286  //
287  // loop over the tracks
288  //
289  if (m_Mode == 1) trackLoopIdT();
290  else trackLoop();
291 
292  if (Debug()) mMiniMcEvent->Print();
293  // Append MC tracks that are not daughters of primary MC vertex.
294  AppendMCDaughterTrack();
295 
296  //
297  // fill the tree and clear all the tracks
298  //
299  if (Debug()) mMiniMcEvent->Print();
300  mMiniMcTree->Fill();
301  mMiniMcEvent->Clear();
302 
303  return kStOk;
304 }
305 
306 /*
307  check if all the association stuff is there
308  */
309 
310 //______________________________________________________________________________
311 Bool_t StMiniMcMaker::initAssociation()
312 {
313  StAssociationMaker* assoc = 0;
314  assoc = (StAssociationMaker*) GetMaker("StAssociationMaker");
315  //
316  // multimaps
317  //
318  if (assoc) {
319  mRcHitMap = assoc->rcTpcHitMap();
320  mRcTrackMap = assoc->rcTrackMap();
321  mMcTrackMap = assoc->mcTrackMap();
322  }
323  return (assoc && mRcHitMap && mRcTrackMap && mMcTrackMap);
324 }
325 
326 /*
327  check if we have a valid vertex
328  */
329 
330 //______________________________________________________________________________
331 Bool_t StMiniMcMaker::initVertex()
332 {
333 
334  if(!mRcEvent->numberOfPrimaryVertices()) {
335  cout << "\tno primary vertex from stevent" << endl;
336  return kFALSE;
337  }
338  if(!mMcEvent->primaryVertex()){
339  cout << "\tno primary vertex from stmcevent" << endl;
340  return kFALSE;
341  }
342  mMainVtx = -1;
343  int maxTks = -1;
344  for (int i=0;i<(int)mRcEvent->numberOfPrimaryVertices();i++) {
345  int numTks = mRcEvent->primaryVertex(i)->numberOfDaughters();
346  if (maxTks > numTks) continue;
347  maxTks = numTks;mMainVtx = i;
348  }
349  if (maxTks <= 0) {
350  cout << "\tEmpty primary vertex from stevent" << endl;
351  return kFALSE;
352  }
353 
354 
355 
356  mRcVertexPos = &mRcEvent->primaryVertex(mMainVtx)->position();
357  mMcVertexPos = &mMcEvent->primaryVertex()->position();
358 
359  if(Debug()){
360  cout
361  << "----------vertex info---------------------\n"
362  << "Position of primary vertex from StEvent: \n"
363  << *mRcVertexPos << endl;
364  cout
365  << "Position of primary vertex from StMcEvent: "<<endl
366  << *mMcVertexPos << endl;
367  cout << "N daughters, StEvent : " << mRcEvent->primaryVertex(mMainVtx)->daughters().size() << endl;
368  cout << "N daughters, StMcEvent : " << mMcEvent->primaryVertex()->daughters().size() << endl;
369  }
370  //
371  // if there was no primary vertex before embedding,
372  // the mc vertex position for each coordinate is equal
373  //
374  return !((std::isnan(mRcVertexPos->x()) || std::isnan(mRcVertexPos->y())));
375 }
376 
377 /*
378 
379  */
380 //______________________________________________________________________________
381 void StMiniMcMaker::trackLoop()
382 {
383  if(Debug()) cout << "##StMiniMcMaker::trackLoop()" << endl;
384 
385  Int_t nMatched(0), nAcceptedRaw(0),nAccepted(0),
386  nMerged(0), nSplit(0), nContam(0), nGhost(0), nMatGlob(0), nContamNew(0),
387  nRcGoodGlobal20(0), nRcGlobal(0), nMcGoodGlobal20(0),
388  nMcNch(0), nMcHminus(0), nMcFtpcWNch(0), nMcFtpcENch(0);
389 
390 
391  RCFOUNDMAP rcFoundMap; // to find split tracks
392  MCFOUNDMAP mcFoundMap; // dont look for a rc match to mc tracks
393  // already flagged as merged or matched
394  //
395  // create the StMiniMcPair class which will hold all
396  // enum Category types.
397  //
398 
399 
400  //
401  // simple loop to associate global tracks
402  // Since the primary track loop has all the bells and
403  // whistles to do the proper accounting of split, merged, contamination
404  // backgrounds, etc. We will do a much simplified version of the
405  // matching loop.
406  // Essentially, I will start only from the monte carlo track container as my seeds,
407  // -I'll apply some simple acceptance cut (10 tpc hits)
408  // -then I will query the association maker for the associated tracks
409  // -I will apply a simple 10 fit points cut to the global tracks
410  // -I will select the one track with the most common hits among those
411  // -enter the matched pair for this global track if there is one.
412  // Note that if there is no associated track, there will be no entry.
413  // The above is all that's needed in the absence of track merging.
414  // To deal with that, I will
415  // -Record which global tracks are already entered
416  // -if a track has been entered, skip the match.
417  // I won't do any more merging accounting, as that is already done for primaries.
418  std::vector<int> enteredGlobalTracks;
419  const StPtrVecMcTrack& allmcTracks = mMcEvent->tracks();
420  cout << "size of mcEvent->tracks() : " << allmcTracks.size() << endl;
421 
422  StMcTrackConstIterator allMcTrkIter = allmcTracks.begin();
423  for ( ; allMcTrkIter != allmcTracks.end(); ++allMcTrkIter) {
424  const StMcTrack* mcGlobTrack = *allMcTrkIter;
425  if(!acceptRaw(mcGlobTrack)) continue; // loose eta cut (4 units, so should include ftpc).
426  if(accept(mcGlobTrack) || mcGlobTrack->ftpcHits().size()>=5) { // 10 tpc hits or 5 ftpc hits
427  if (Debug()>1) {
428  cout << "accepted mc global track, key " << mcGlobTrack->key() << endl;
429  }
430  // Ok, track is accepted, query the map for its reco tracks.
431 
432  StTrackPairInfo* candTrackPair = findBestMatchedGlobal(mcGlobTrack);
433 
434  if (candTrackPair) {
435  // ok, found a match! Enter into the array and store the glob id
436  // to only enter a glob track once
437  const StGlobalTrack* glTrack = candTrackPair->partnerTrack();
438 
439  if (Debug()>1) {
440  cout << "global match " << endl;
441  cout << "mc, rc pt " << mcGlobTrack->pt() << ", " << glTrack->geometry()->momentum().perp() << endl;
442  }
443 
444  if (find(enteredGlobalTracks.begin(),enteredGlobalTracks.end(),glTrack->key())!=enteredGlobalTracks.end()) continue; //if it's already matched, skip it.
445  StMiniMcPair* miniMcPair = new StMiniMcPair;
446  Int_t commonHits =
447  candTrackPair->commonTpcHits()%100+
448  ((candTrackPair->commonSvtHits()%10)*100)+
449  ((candTrackPair->commonSsdHits()%10)*1000);
450  fillTrackPairInfo(miniMcPair, mcGlobTrack,
451  0, glTrack,
452  commonHits, mRcTrackMap->count(glTrack),
453  mMcTrackMap->count(mcGlobTrack), 0,
454  kTRUE);
455  mMiniMcEvent->addTrackPair(miniMcPair,MATGLOB);
456  delete miniMcPair;
457  enteredGlobalTracks.push_back(glTrack->key()); // store the keys of the tracks we've matched.
458  nMatGlob++;
459 
460  }
461  }// mc hits condition
462  }// end of global track match loop
463 
464  // primary track begins here
465  //
466  // loop over mc tracks.
467  //
468 
469  const StPtrVecMcTrack& mcTracks = mMcEvent->primaryVertex()->daughters();
470 
471  cout << "size of MC primary tracks : " << mcTracks.size() << endl;
472 
473  StMcTrackConstIterator mcTrkIter = mcTracks.begin();
474  for( ; mcTrkIter != mcTracks.end(); mcTrkIter++){
475  const StMcTrack* mcTrack = *mcTrkIter;
476 
477  // DEBUG
478  if(!mcTrack) { cout << "No mc track? " << endl; continue; }
479 
480  //if(!acceptPt(mcTrack)) continue; // pt cut?
481  if(!acceptRaw(mcTrack)) continue; // loose eta cuts, etc
482 
483  nAcceptedRaw++;
484 
485  if(fabs(mcTrack->pseudoRapidity())<.5 && isPrimaryTrack(mcTrack) && mcTrack->particleDefinition()){
486  if(mcTrack->particleDefinition()->charge()!=0) nMcNch++;
487  if(mcTrack->particleDefinition()->charge()<0) nMcHminus++;
488  }
489  if(mcTrack->particleDefinition() && mcTrack->particleDefinition()->charge()!=0 && isPrimaryTrack(mcTrack) ) {
490  if(mcTrack->pseudoRapidity()<-2.8 && mcTrack->pseudoRapidity()>-3.8) nMcFtpcENch++;
491  if(mcTrack->pseudoRapidity()>2.8 && mcTrack->pseudoRapidity()<3.8) nMcFtpcWNch++;
492  }
493 
494  if(acceptGood20(mcTrack)) nMcGoodGlobal20++;
495 
496  Int_t nAssocGl = mMcTrackMap->count(mcTrack);
497  Int_t nAssocPr = 0; // value maybe reset below.
498 
499  //
500  // minimum requirement to accept the mc track and search a rc match
501  //
502  if(accept(mcTrack)){ //10 tpc hits
503 
504  nAccepted++;
505 
506  // check if already included as merged or matched
507  //
508  if(!mcFoundMap.count(mcTrack->key())){
509 
510  //
511  // follow manuel's example.
512  // find the best matched associated rc track only.
513  //
514 
515  PAIRVEC candPair = findMatchedRc(mcTrack);
516  nAssocPr = candPair.size(); // # of primaries matched to mcTrack
517 
518  if(candPair.size()>0){ // match to at least one primary track
519 
520  // find the best matched rc track
521 
522  PAIRVEC::iterator iterBestMatchPair
523  = max_element(candPair.begin(),candPair.end(),pairCmp);
524 
525  //*** DEBUG
526  /*
527  if(candPair.size()>1){
528  cout << "best:" << endl
529  << "\t" << (*iterBestMatchPair)->commonTpcHits() << endl;
530  cout << "all: " << endl;
531  for(unsigned int j=0; j<candPair.size(); j++){
532  cout << "\t" << candPair[j]->commonTpcHits() << endl;
533  }
534  }
535  */
536  //*** DEBUGEND
537 
538 
539  const StGlobalTrack* glTrack = (*iterBestMatchPair)->partnerTrack();
540  const StPrimaryTrack* prTrack = isPrimaryTrack(glTrack);
541 
542  //
543  // 02/15/01
544  // safer idea. given an mc track, find the pr track
545  // that's best matched. but then find and sort all the mc tracks
546  // that are best matched to the same pr track.
547  // the best of these best matched mc tracks will be stored
548  // as 'matched', while the other will be stored as merged.
549  // the user can then apply a tighter(or looser) common hits
550  // cut downstream to reduce these # of merged.
551  // probably much ado about nothing since the number of
552  // merged is very small.
553  //
554 
555  PAIRVEC mcMergedPair; // actually, merged or matched pairs
556 
557  pair<rcTrackMapIter,rcTrackMapIter> rcBounds
558  = mRcTrackMap->equal_range(glTrack);
559  rcTrackMapIter rcMapIter = rcBounds.first;
560 
561  // # of mc tracks matched to this rc track
562  //
563  UInt_t nAssocMc = mRcTrackMap->count(glTrack);
564  std::vector<UInt_t> nAssocGlVec; // # of rc globals matched to a merged mc cand
565  std::vector<UInt_t> nAssocPrVec; // # of rc primaries matched to a merged mc
566  //
567  // loop over the mc tracks associated with this rc track
568  //
569 
570  for( ; rcMapIter != rcBounds.second; rcMapIter++){
571  StTrackPairInfo* assocPair = (*rcMapIter).second;
572  const StMcTrack* mcCandTrack = assocPair->partnerMcTrack();
573 
574  // no pt cut
575  // but must lie with the loose cut (acceptRaw)
576  // and mc hits cut (accept)
577  if(!acceptRaw(mcCandTrack) || !accept(mcCandTrack)) continue;
578 
579  //
580  // loop over the rc tracks matched with this cand mc track.
581  // is this mc track best matched with the original rc track?
582  //
583  PAIRVEC rcCandPair = findMatchedRc(mcCandTrack);
584 
585  PAIRVEC::iterator iterRcBestPair =
586  max_element(rcCandPair.begin(),rcCandPair.end(),pairCmp);
587 
588  // yes. this is a merged track candidate.
589  //
590  if((*iterRcBestPair)->partnerTrack() == glTrack){
591 
592  nAssocGlVec.push_back(mMcTrackMap->count(mcCandTrack));
593  nAssocPrVec.push_back(rcCandPair.size());
594 
595  mcMergedPair.push_back(assocPair);
596  }
597  }
598 
599  //
600  // now we have all the possible mc merged tracks to this rc track.
601  // (always sort. probably not cpu intensive to sort one element)
602  //
603  sort(mcMergedPair.begin(),mcMergedPair.end(),sortCmp);
604 
605  //
606  // ok, save all primary mc tracks.
607  // best 'best' matched is saved as a matched pair.
608 
609  if(Debug()==2 && mcMergedPair.size()>1) {
610  cout << "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%" << endl;
611  cout << "MERGED" << endl;
612  }
613 
614  Bool_t foundBest = kFALSE;
615  Bool_t isBestContam = kFALSE; // the best 'best matched' mc track
616  // is not a primary mc
617  for(unsigned int i=0; i<mcMergedPair.size(); i++){
618  const StMcTrack* mergedMcTrack = (mcMergedPair[i])->partnerMcTrack();
619  UInt_t mergedCommonHits = (mcMergedPair[i])->commonTpcHits();
620 
621  if(Debug()==2 && mcMergedPair.size()>1) {
622  TString hello = (foundBest) ? "yes" : "no";
623  cout << "-----------" << " foundBest? " << hello.Data() << endl;
624  checkMerged(mergedMcTrack, mergedCommonHits,prTrack);
625  }
626 
627  // only primary mc
628  if(!isPrimaryTrack(mergedMcTrack)){
629  if(i==0) isBestContam = kTRUE;
630  continue;
631  }
632 
633  if(!foundBest){
634  // 02/02/02 rc pt cut
635  if(acceptPt(glTrack) || acceptPt(prTrack)){
636  StMiniMcPair* miniMcPair = new StMiniMcPair;
637  Int_t commonHits = mergedCommonHits%100+
638  ((mcMergedPair[i]->commonSvtHits()%10)*100)+
639  ((mcMergedPair[i]->commonSsdHits()%10)*1000);
640  fillTrackPairInfo(miniMcPair, mergedMcTrack,
641  prTrack, glTrack,
642  commonHits, nAssocMc,
643  nAssocGlVec[i], nAssocPrVec[i],
644  isBestContam);
645  mMiniMcEvent->addTrackPair(miniMcPair,MATCHED);
646  delete miniMcPair;
647  }
648  rcFoundMap[prTrack->key()]=1; // the value is meaningless
649  nMatched++;
650  foundBest = kTRUE;
651  }
652  else{
653  // 02/02/02 rc pt cut
654  if(acceptPt(glTrack) || acceptPt(prTrack)){
655  StMiniMcPair* miniMcPair = new StMiniMcPair;
656  Int_t commonHits =
657  mergedCommonHits+
658  ((mcMergedPair[i]->commonSvtHits()%10)*100)+
659  ((mcMergedPair[i]->commonSsdHits()%10)*1000);
660  fillTrackPairInfo(miniMcPair,mergedMcTrack,prTrack,glTrack,
661  commonHits, nAssocMc,nAssocGlVec[i],
662  nAssocPrVec[i]);
663  mMiniMcEvent->addTrackPair(miniMcPair,MERGED);
664  delete miniMcPair;
665  }
666  if(Debug()==2 && acceptDebug(mergedMcTrack))
667  cout << "YES! satisfies cuts" << endl;
668 
669  nMerged++;
670 
671  }
672 
673  //
674  // flag this mc track so we dont process it again.
675  //
676  mcFoundMap[mergedMcTrack->key()]=1;
677 
678  } // 'merged' pair loop
679  if(Debug()==2 && mcMergedPair.size()>1)
680  cout << "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%" << endl;
681  } // found a matched rc track
682  } // not merged
683  } // accept mc track
684 
685  //
686  // fill the RAW mc track info here, since we need the number of rc
687  // primaries associated with this mc track
688  //
689  // 02/25/02
690  // accept all mc primary tracks within some eta window
691  if(acceptRaw(mcTrack)){
692  StTinyMcTrack tinyMcTrack;
693  fillMcTrackInfo(&tinyMcTrack,mcTrack,nAssocGl,nAssocPr);
694  mMiniMcEvent->addMcTrack(&tinyMcTrack);
695  }
696  } // mc track iter
697 
698  //nSplit(0), nContam(0), nGhost(0), nContamNew(0);
699 
700  //
701  // need to loop over the rc primary tracks to get
702  // the event centrality values, etc.
703  // also look for 'split tracks','ghost tracks'...
704  //
705 
706  Int_t nUncorrectedGlobals(0);
707  StSPtrVecTrackNode& trackNode = mRcEvent->trackNodes();
708  int nTracks = trackNode.size();
709  for (int i=0; i < nTracks; i++) {
710  StTrackNode *node = trackNode[i];
711  if (!node) continue;
712  StTrack *glTrack = node->track(global);
713  if (! glTrack) continue;
714  if(acceptGlobals(glTrack)) nUncorrectedGlobals++;
715  }
716 
717  Int_t nGoodTrackEta(0), nFtpcWUncorrected(0), nFtpcEUncorrected(0);
718 
719  const StSPtrVecPrimaryTrack& prTracks =
720  mRcEvent->primaryVertex(mMainVtx)->daughters();
721 
722  for(UInt_t i=0; i<prTracks.size(); i++){
723  StPrimaryTrack* prTrack = prTracks[i];
724 
725  //
726  // check for positive flag
727  //
728  if(!ok(prTrack)) continue;
729 
730  const StGlobalTrack* glTrack
731  = static_cast<const StGlobalTrack*>(prTrack->node()->track(global));
732 
733 
734  //
735  // centrality
736  //
737  if(acceptCentrality(prTrack)) nGoodTrackEta++;
738  nRcGlobal++; //These are identical and ok(prTrack) checks the flag
739  if(acceptGood20(glTrack)) nRcGoodGlobal20++;
740 
741  //
742  // uncorrected negative primaries
743  // 02/25/02 no longer used. uses manuels function in fillEventInfo
744  if(acceptFTPC(prTrack) && prTrack->geometry()->momentum().pseudoRapidity() < 0. ) nFtpcEUncorrected++;
745  if(acceptFTPC(prTrack) && prTrack->geometry()->momentum().pseudoRapidity() > 0. ) nFtpcWUncorrected++;
746 
747  //
748  // rc pt cut ?
749  //
750  if(!acceptPt(glTrack) && !acceptPt(prTrack)) continue;
751 
752  //
753  // minimum cut
754  //
755  if(!accept(glTrack) || !accept(prTrack)) continue;
756 
757  // tracks we'll put in the tree.
758  mNRc++; // for printing in Finish()
759 
760  UInt_t nAssocGl=0, nAssocPr=0;
761  UInt_t nAssocMc = mRcTrackMap->count(glTrack);
762 
763  //
764  // tracks that were matched to a mc primary but never best matched..
765  // split?
766  //
767  if(mRcTrackMap->count(glTrack) && !rcFoundMap.count(prTrack->key())){
768 
769  // make sure the match was best matched to a mc primary
770  //
771 
772  pair<rcTrackMapIter,rcTrackMapIter> rcBounds
773  = mRcTrackMap->equal_range(glTrack);
774 
775  PAIRVEC candPair;
776 
777  rcTrackMapIter rcMapIter = rcBounds.first;
778 
779  for( ; rcMapIter != rcBounds.second; rcMapIter++){
780  StTrackPairInfo* assocPair = (*rcMapIter).second;
781  //const StMcTrack* mcCandTrack = assocPair->partnerMcTrack();
782 
783  //
784  // new bug. make sure the mc track was accepted!
785  //
786  // no pt cut
787  //** 01/28/02
788  // now there's absolutely no cut on the mc track.
789  // if(!acceptRaw(mcCandTrack) || !accept(mcCandTrack)) continue;
790 
791  candPair.push_back(assocPair); // for contamination
792 
793  }
794 
795  //** 01/28/02 - now this should always be true.
796  if(candPair.size()){ // found a matching mc track that passes the cuts
797 
798  PAIRVEC::iterator iterBestMatchPair
799  = max_element(candPair.begin(),candPair.end(),pairCmp);
800 
801  const StMcTrack* mcTrack = (*iterBestMatchPair)->partnerMcTrack();
802  UInt_t commonHits = (*iterBestMatchPair)->commonTpcHits();
803 
804  nAssocGl = mMcTrackMap->count(mcTrack);
805 
806  //
807  // how many rc primaries does this mc track match with?
808  //
809  PAIRVEC matchPair = findMatchedRc(mcTrack);
810  nAssocPr = matchPair.size();
811 
812  //
813  // best matched to a mc primary ?
814  // yes, it's a split track
815  //
816  if(isPrimaryTrack(mcTrack)){
817  StMiniMcPair* miniMcPair = new StMiniMcPair;
818  Int_t cHits =
819  commonHits%100+
820  (((*iterBestMatchPair)->commonSvtHits()%10)*100)+
821  (((*iterBestMatchPair)->commonSsdHits()%10)*1000);
822  fillTrackPairInfo(miniMcPair,mcTrack,prTrack,glTrack,
823  cHits, nAssocMc, nAssocGl, nAssocPr);
824  mMiniMcEvent->addTrackPair(miniMcPair,SPLIT);
825  delete miniMcPair;
826  nSplit++;
827 
828  //
829  // reality check.
830  //
831  if(Debug()==2) checkSplit(mcTrack,glTrack,commonHits);
832  }
833  else{ // no, it's best matched to a non primary, contamination
834  StContamPair* contamPair = new StContamPair;
835  Int_t cHits =
836  commonHits%100+
837  (((*iterBestMatchPair)->commonSvtHits()%10)*100)+
838  (((*iterBestMatchPair)->commonSsdHits()%10)*1000);
839  fillTrackPairInfo(contamPair,mcTrack,
840  prTrack,glTrack,cHits,
841  nAssocMc,nAssocGl,nAssocPr);
842  mMiniMcEvent->addContamPair(contamPair);
843  delete contamPair;
844  if(Debug()>=2) checkContam(mcTrack,glTrack,commonHits);
845  nContam++;
846  }
847  } // cand.size()?
848  else{
849  cout << "ERROR- no match for rc track when determining split/contaminations? " << endl;
850  exit(1);
851  }
852 
853  } // split?
854  //
855  // ghosts?
856  //
857  else if(mGhost && mRcTrackMap->count(glTrack)==0){
858  StMiniMcPair* miniMcPair = new StMiniMcPair;
859  fillRcTrackInfo(miniMcPair,
860  prTrack,glTrack,nAssocMc);
861  mMiniMcEvent->addTrackPair(miniMcPair,GHOST);
862  delete miniMcPair;
863  nGhost++;
864  if(Debug()) {
865  cout << "#############" << endl;
866  cout << "GHOST!" << endl;
867  cout << "pr pt: " << prTrack->geometry()->momentum().perp() << endl
868  << "TPC fit hits : " << glTrack->fitTraits().numberOfFitPoints(kTpcId)
869  << "SVT fit hits : " << glTrack->fitTraits().numberOfFitPoints(kSvtId)
870  << endl;
871  }
872  }
873  }
874  //
875  // fill all the event info
876  //
877 
878  fillEventInfo(nGoodTrackEta,nRcGlobal,nRcGoodGlobal20,
879  nAcceptedRaw,
880  nMcGoodGlobal20,nMcNch,nMcHminus,
881  nMcFtpcENch, nMcFtpcWNch,nFtpcEUncorrected,nFtpcWUncorrected,
882  nUncorrectedGlobals);
883 
884  // delete the most probable pid functor
885  //
886  // delete mPidAlgo;
887 
888  cout << "\tall rc tracks: " << prTracks.size() << endl;
889  cout << "\tn good eta : " << nGoodTrackEta << endl;
890  cout << "\tcentrality : " << mMiniMcEvent->centrality() << endl;
891  cout << "\tuncorrected : " << mMiniMcEvent->nUncorrectedPrimaries() << endl;
892  cout << "\tall mc tracks: " << mcTracks.size() << endl;
893  cout << "\taccepted raw : " << nAcceptedRaw << endl;
894  cout << "\taccepted mc : " << nAccepted << endl;
895  cout << "\tmatched rc gl: " << nMatGlob << endl;
896  cout << "\tmatched rc pr: " << nMatched << endl;
897  cout << "\tmerged rc : " << nMerged << endl;
898  cout << "\tsplit rc : " << nSplit << endl;
899  if(mGhost) {
900  cout << "\tghost rc : " << nGhost << endl;
901  cout << "\tcontam rc : " << nContam << endl;
902  cout << "\tcontam new : " << nContamNew << endl;
903  }
904 
905  // counters
906  mNSplit += nSplit; mNGhost+= nGhost; mNContam += nContam;
907  mNMatched += nMatched;
908  mNMatGlob += nMatGlob;
909 }
910 //________________________________________________________________________________
911 struct MyHolder_t { const StTrack *gl;const StTrack *pr;int hits; float qa;};
912 void StMiniMcMaker::trackLoopIdT() // match with IdTruth
913 {
914 
915 
916  typedef std::map< int,StTinyMcTrack*> McTinyMap_t;
917  typedef McTinyMap_t::iterator McTinyMapIter_t;
918 
919 // struct MyHolder_t { public: const StTrack *gl;const StTrack *pr;int hits; float qa;};
920 
921  typedef std::multimap< int,MyHolder_t> StTrackMap_t;
922  typedef std::pair<int,MyHolder_t> StTrackPair_t;
923  typedef StTrackMap_t::const_iterator StTrackMapIter_t;
924 
925  Int_t nAcceptedRaw(0),nAccepted(0),
926  nRcGoodGlobal20(0), nRcGlobal(0), nMcGoodGlobal20(0),
927  nMcNch(0), nMcHminus(0), nMcFtpcWNch(0), nMcFtpcENch(0);
928  Int_t nGoodTrackEta(0), nFtpcWUncorrected(0), nFtpcEUncorrected(0);
929  Int_t nUncorrectedGlobals(0);
930 
931  const StVertex *rcVertex = mRcEvent->primaryVertex(mMainVtx);
932 
933  const StPtrVecMcTrack& mcTracks = mMcEvent->tracks();
934  int NMcTracks = mcTracks.size();
935  cout << "size of mcEvent->tracks() : " << NMcTracks << endl;
936 
937 // Fill global and primary StTrackMap's
938  StTrackMap_t glMap,prMap;
939  StSPtrVecTrackNode& trackNode = mRcEvent->trackNodes();
940  int nTracks = trackNode.size(),myKey;
941  for (int i=0; i < nTracks; i++) {
942  StTrackNode *node = trackNode[i];
943  if (!node) continue;
944  StTrack *glTrack = node->track(global);
945  if (! glTrack) continue;
946  if(acceptGlobals(glTrack)) nUncorrectedGlobals++;
947 
948  nRcGlobal++;
949  if(acceptGood20(glTrack)) nRcGoodGlobal20++;
950 
951  StPrimaryTrack *prTrack = (StPrimaryTrack*)node->track(primary);
952  if (prTrack && prTrack->vertex()!=rcVertex) prTrack=0;
953  MyHolder_t h; h.gl = glTrack;h.pr = prTrack;
954  dominatTkInfo(glTrack,myKey ,h.hits,h.qa);
955  glMap.insert(StTrackPair_t(myKey,h));
956  if (!prTrack) continue;
957  prMap.insert(StTrackPair_t(myKey,h));
958 
959  if(acceptCentrality(prTrack)) nGoodTrackEta++;
960  if(acceptFTPC(prTrack) && prTrack->geometry()->momentum().pseudoRapidity() < 0. ) nFtpcEUncorrected++;
961  if(acceptFTPC(prTrack) && prTrack->geometry()->momentum().pseudoRapidity() > 0. ) nFtpcWUncorrected++;
962  }
963 
964 // Fill Tiny Mc tracks
965  McTinyMap_t mcTinyMap ;
966  for (int k = 0; k < NMcTracks; k++) {
967  const StMcTrack* mcTrack = mcTracks[k];
968  if (! mcTrack->geantId()) continue;
969  int tkKey = mcTrack->key(); if (!tkKey) continue;
970  if(!acceptRaw(mcTrack)) continue; // loose eta cuts, etc
971  nAcceptedRaw++;
972  if(!accept(mcTrack)) continue;
973  nAccepted++;
974  if(acceptGood20(mcTrack)) nMcGoodGlobal20++;
975  if(fabs(mcTrack->pseudoRapidity())<.5 && isPrimaryTrack(mcTrack) && mcTrack->particleDefinition()){
976  if(mcTrack->particleDefinition()->charge()!=0) nMcNch++;
977  if(mcTrack->particleDefinition()->charge() <0) nMcHminus++;
978  }
979  if(mcTrack->particleDefinition() && mcTrack->particleDefinition()->charge()!=0 && isPrimaryTrack(mcTrack) ) {
980  if(mcTrack->pseudoRapidity()<-2.8 && mcTrack->pseudoRapidity()>-3.8) nMcFtpcENch++;
981  if(mcTrack->pseudoRapidity()> 2.8 && mcTrack->pseudoRapidity()< 3.8) nMcFtpcWNch++;
982  }
983  int nAssocGl = glMap.count(tkKey);
984  int nAssocPr = prMap.count(tkKey);
985  StTinyMcTrack *tinyMcTrack = mMiniMcEvent->addMcTrack(0);
986  fillMcTrackInfo(tinyMcTrack,mcTrack,nAssocGl,nAssocPr);
987  mcTinyMap[tkKey] = tinyMcTrack;
988  tinyMcTrack->ResetBit(1);
989  }
990 
991 // Fill miniMcPair's by matched info
992 
993  for (StTrackMapIter_t it =glMap.begin(); it!= glMap.end();++it) {
994  int tkKey = (*it).first;
995  const MyHolder_t &h = (*it).second;
996  StMiniMcPair *miniMcPair = mMiniMcEvent->addTrackPair(0,MATCHED);
997  StTinyMcTrack *tinyMcTrack = mcTinyMap[tkKey];
998  if (tinyMcTrack) {
999  *(StTinyMcTrack*)miniMcPair = *tinyMcTrack;
1000  tinyMcTrack->SetBit(1); //Mark this as used
1001  }
1002  fillRcTrackInfo(miniMcPair,h.pr,h.gl,0);
1003  miniMcPair->setDominatrack(tkKey);
1004  miniMcPair->setDominCommonHit(h.hits);
1005  miniMcPair->setAvgQuality(h.qa);
1006  miniMcPair->setNCommonHit(h.hits); //??VP
1007  }
1008 
1009 //
1010 // Fill miniMcPair's by UNmatched info
1011  for( McTinyMapIter_t it =mcTinyMap.begin();it !=mcTinyMap.end();++it)
1012  {
1013  StTinyMcTrack *tinyMcTrack = (*it).second;
1014  if (!tinyMcTrack) continue;
1015  if (tinyMcTrack->TestBit(1)) continue;
1016  StMiniMcPair *miniMcPair = mMiniMcEvent->addTrackPair(0,MATCHED);
1017  *(StTinyMcTrack*)miniMcPair = *tinyMcTrack;
1018  }
1019 
1020 
1021 
1022 
1023  fillEventInfo(nGoodTrackEta,nRcGlobal,nRcGoodGlobal20,
1024  nAcceptedRaw,nMcGoodGlobal20,nMcNch ,nMcHminus,
1025  nMcFtpcENch ,nMcFtpcWNch ,nFtpcEUncorrected,nFtpcWUncorrected,
1026  nUncorrectedGlobals);
1027 }
1028 //______________________________________________________________________________
1029 void StMiniMcMaker::buildEmcIndexArray()
1030 {
1031  if (Debug()>1) cout << "StMiniMcMaker::buildEmcIndexArray" << endl;
1032  // Put all the EMC hits into an array that can then be
1033  // indexed by SoftId, for use when storing tower info for
1034  // tracks.
1035 
1036  // zero all pointers in the array
1037  // Note: vector::clear() is not what should be used
1038  // clear() removes all the elements from the vector, i.e. the
1039  // size() after a clear() will be zero.
1040  for (int softIdNum=1; softIdNum<4801; ++softIdNum) mEmcIndex[softIdNum]=0;
1041 
1042  StEmcGeom* emcGeom = StEmcGeom::getEmcGeom(1);
1043  if (! mRcEvent->emcCollection()) return;
1044  StEmcDetector* bemcDet = mRcEvent->emcCollection()->detector(kBarrelEmcTowerId);
1045  if (! bemcDet) return;
1046  if (Debug()>1) {
1047  cout << "emcGeom " << emcGeom << endl;
1048  cout << "bemcDet " << bemcDet << endl;
1049  }
1050  // Loop over modules. Note, the modules are returned by module number
1051  // i.e. in the range 1-120 for BEMC.
1052  for (size_t iMod=1; iMod<=bemcDet->numberOfModules(); ++iMod) {
1053  StSPtrVecEmcRawHit& modHits = bemcDet->module(iMod)->hits();
1054  if (Debug()>1) {
1055  cout << "Module " << iMod << endl;
1056  cout << "Hits in Module " << modHits.size() << endl;
1057  }
1058  for (size_t iHit=0; iHit<modHits.size(); ++iHit) {
1059  StEmcRawHit* rawHit = modHits[iHit];
1060  unsigned int softId = rawHit->softId(1); // 1 is "bemc" (see StEmcGeom::getDetNumFromName)
1061  if (Debug()>2) {
1062  cout << "iHit " << iHit << endl;
1063  cout << "soft Id " << softId << endl;
1064  cout << "ene, mod, eta, sub " << rawHit->energy()
1065  << ", " << rawHit->module()
1066  << ", " << rawHit->eta()
1067  << ", " << rawHit->sub()
1068  << endl;
1069  }
1070  // Make sure the softId is valid.
1071  // Note: checkId returns 0 (false) when the softId is between 1 and mNRaw
1072  // so need an extra "!" to use inside the if statement.
1073  if (!emcGeom->checkId(softId)) {
1074  mEmcIndex[softId] = rawHit;
1075  } //check for valid ID
1076  } // loop over hits
1077  } // loop over modules
1078  return;
1079 }
1080 
1081 //________________________________________________________________________________
1082 /*
1083  Create the output root file and the TTree
1084 */
1085 
1086 Int_t StMiniMcMaker::openFile()
1087 {
1088  cout << "###StMiniMcMaker::openFile()" << endl;
1089 
1090  //
1091  // for the output root file,
1092  // - replace geant.root with minimc.root
1093  // in case there is an event.root, or any other blabla.root, we will
1094  // replace the blabla.root with minimc.root. This takes some funky TString
1095  // manipulation, and will only work if the filename is of the form somefilestring.blablabla.root
1096  // it will change it to somefilestring.minimc.root
1097  // - for the case of running in bfc, the Infilename will have only the extension .root
1098  // so we can keep this working in both cases by checking that the first and the last "." are
1099  // not the same before attempting to remove what's in between them.
1100  cout << "Infilename = " << mInFileName << endl;
1101  if (mInFileName == "") {
1102  TFile *tfile = GetTFile();
1103  assert(tfile);
1104  tfile->cd();
1105  } else {
1106  TString outFileName(mInFileName);
1107  //outFileName.ReplaceAll("geant.root","minimc.root");
1108  short indx1 = outFileName.First('.');
1109  short indx2 = outFileName.Last('.');
1110  if (indx1!=indx2) outFileName.Remove(indx1+1,(indx2-indx1));
1111  outFileName.Insert(indx1+1,"minimc.");
1112  outFileName.Prepend(mOutDir);
1113 
1114  closeFile();
1115  cout << "Opening File " << outFileName << endl;
1116  mMiniMcDST = new TFile(outFileName.Data(),"RECREATE");
1117 
1118  if(!mMiniMcDST || !mMiniMcDST->IsOpen()){
1119  gMessMgr->Error() << "Cannot create = " << outFileName << endm;
1120  return kStErr;
1121  }
1122  // mMiniMcDST->SetFormat(1);
1123  cout << "\toutfile = " << outFileName << endl;
1124  mOutFileName = outFileName;
1125  //
1126  // top level tree
1127  //
1128  }
1129  if(Debug()) cout << "##Creating the top level tree..." << endl;
1130  mMiniMcTree = new TTree("StMiniMcTree","StMiniMcTree");
1131 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,01,05)
1132  mMiniMcTree->SetBranchStyle(0);
1133 #endif
1134  Int_t bufSZ = 64000; // make this bigger?
1135  mMiniMcTree->Branch("StMiniMcEvent","StMiniMcEvent",&mMiniMcEvent, bufSZ,1);
1136  mMiniMcTree->SetAutoSave(10000000); // 10 MB
1137 
1138  cout << "##...done" << endl;
1139 
1140  return kStOk;
1141 }
1142 
1143 //______________________________________________________________________________
1144 // close file and write.
1145 Int_t StMiniMcMaker::closeFile()
1146 {
1147  cout << "###StMiniMcMaker::closeFile()" << endl;
1148  cout << "\tWriting " << mOutFileName << endl;
1149 
1150  if(mMiniMcDST && mMiniMcDST->IsOpen()){
1151  mMiniMcDST->cd();
1152  mMiniMcTree->Write();
1153  mMiniMcDST->Write();
1154  mMiniMcDST->Close();
1155  delete mMiniMcDST; mMiniMcDST=0;
1156  }
1157  mMiniMcTree=0; //It is deleted by TFile::Close
1158  cout << "\t...done\n";
1159 
1160  return kStOk;
1161 }
1162 
1163 //______________________________________________________________________________
1164 void StMiniMcMaker::fillEventInfo(Int_t nGoodTrackEta, Int_t nRcGlobal, Int_t nRcGoodGlobal20,
1165  Int_t nMcGlobal, Int_t nMcGoodGlobal20,
1166  Int_t nMcNch, Int_t nMcHminus, Int_t nMcFtpcENch, Int_t nMcFtpcWNch, Int_t nFtpcEUncorrected,
1167  Int_t nFtpcWUncorrected, Int_t nUncorrectedGlobals)
1168 {
1169  mMiniMcEvent->setEventId((Int_t) mRcEvent->id());
1170  mMiniMcEvent->setRunId((Int_t) mRcEvent->runId());
1171  mMiniMcEvent->setOriginMult((Int_t)mRcEvent->primaryVertex(mMainVtx)->numberOfDaughters());
1172  mMiniMcEvent->setCentralMult(nGoodTrackEta);
1173 
1174  mMiniMcEvent->setImpact (mMcEvent->impactParameter() );
1175  mMiniMcEvent->setImpactPhi (mMcEvent->phiReactionPlane() );
1176  mMiniMcEvent->setTimeOffset (mMcEvent->triggerTimeOffset());
1177 
1178  mMiniMcEvent->setNMcNch(nMcNch);
1179  mMiniMcEvent->setNMcFtpcWNch(nMcFtpcWNch);
1180  mMiniMcEvent->setNMcFtpcENch(nMcFtpcENch);
1181  mMiniMcEvent->setNMcHminus(nMcHminus);
1182 
1183  mMiniMcEvent->setNMcGlobal(nMcGlobal); // from nAcceptedRaw, no point in doing two variables for the same number
1184  mMiniMcEvent->setNMcGoodGlobal20(nMcGoodGlobal20);
1185 
1186  mMiniMcEvent->setNRcGlobal(nRcGlobal);
1187  mMiniMcEvent->setNRcGoodGlobal20(nRcGoodGlobal20);
1188 
1189  mMiniMcEvent->setNUncorrectedNegativePrimaries(uncorrectedNumberOfNegativePrimaries(*mRcEvent));
1190  mMiniMcEvent->setNUncorrectedPrimaries(uncorrectedNumberOfPrimaries(*mRcEvent));
1191 
1192  mMiniMcEvent->setNUncorrectedGlobals(nUncorrectedGlobals);
1193 
1194  mMiniMcEvent->setNFtpcWUncorrectedPrimaries(nFtpcWUncorrected);
1195  mMiniMcEvent->setNFtpcEUncorrectedPrimaries(nFtpcEUncorrected);
1196 
1197  mMiniMcEvent->setCentrality(getIndex((size_t) mMiniMcEvent->nUncorrectedPrimaries()));
1198  mMiniMcEvent->setMcMult(mMcEvent->numberOfPrimaryTracks());
1199 
1200  mMiniMcEvent->setVertexX(mRcVertexPos->x());
1201  mMiniMcEvent->setVertexY(mRcVertexPos->y());
1202  mMiniMcEvent->setVertexZ(mRcVertexPos->z());
1203  StMatrixF C(mRcEvent->primaryVertex(mMainVtx)->covariantMatrix());
1204  Float_t cov[6] = {C(1,1),C(2,1),C(2,2),C(3,1),C(3,2),C(3,3)};
1205  mMiniMcEvent->setVertexCovMatrix(cov);
1206  mMiniMcEvent->setMcVertexX(mMcVertexPos->x());
1207  mMiniMcEvent->setMcVertexY(mMcVertexPos->y());
1208  mMiniMcEvent->setMcVertexZ(mMcVertexPos->z());
1209 
1210  if (mRcEvent->runInfo()) {
1211  mMiniMcEvent->setMagField(static_cast<Float_t>(mRcEvent->runInfo()->magneticField()));
1212  mMiniMcEvent->setBackgroundRate(mRcEvent->runInfo()->backgroundRate());
1213  mMiniMcEvent->setCenterOfMassEnergy(mRcEvent->runInfo()->centerOfMassEnergy());
1214  mMiniMcEvent->setBeamMassNumberEast(mRcEvent->runInfo()->beamMassNumber(east));
1215  mMiniMcEvent->setBeamMassNumberWest(mRcEvent->runInfo()->beamMassNumber(west));
1216  }
1217 
1218  Float_t ctb = -1., zdce = -1, zdcw = -1;
1219 
1220  StTriggerDetectorCollection *triggers
1221  = mRcEvent->triggerDetectorCollection();
1222  if (triggers) {
1223  StCtbTriggerDetector &CTB = triggers->ctb();
1224  StZdcTriggerDetector &ZDC = triggers->zdc();
1225  // get CTB
1226  for (UInt_t slat=0; slat<CTB.numberOfSlats(); slat++) {
1227  for (UInt_t tray=0; tray<CTB.numberOfTrays();tray++) {
1228  ctb += CTB.mips(tray,slat,0);
1229  }
1230  }
1231  //get ZDCe and ZDCw
1232  zdce = ZDC.adcSum(east);
1233  zdcw = ZDC.adcSum(west);
1234  }
1235 
1236  mMiniMcEvent->setCtb(ctb);
1237  mMiniMcEvent->setZdcE(zdce);
1238  mMiniMcEvent->setZdcW(zdcw);
1239 
1240 }
1241 
1242 //______________________________________________________________________________
1243 void StMiniMcMaker::fillTrackPairInfo( StMiniMcPair* miniMcPair,
1244  const StMcTrack* mcTrack,
1245  const StTrack* prTrack,
1246  const StTrack* glTrack,
1247  Int_t commonHits,
1248  Int_t nAssocMc, Int_t nAssocGl,
1249  Int_t nAssocPr, Bool_t isBestContam)
1250 {
1251 
1252  if(mcTrack) fillMcTrackInfo(miniMcPair,mcTrack,nAssocGl,nAssocPr);
1253 
1254  if(prTrack || glTrack) fillRcTrackInfo(miniMcPair,prTrack,glTrack,nAssocMc);
1255 
1256  // common association info
1257  miniMcPair->setNCommonHit(commonHits);
1258  miniMcPair->setIsBestContam(isBestContam);
1259  //
1260  int aeonFlux(-999);
1261  short aeonFluxHits(0);
1262  float aeonFluxQuality(-999); // MCBS: name suggested by Jerome.
1263  if (prTrack) dominatrackInfo(prTrack,aeonFlux,aeonFluxHits,aeonFluxQuality);
1264  else if (glTrack) dominatrackInfo(glTrack,aeonFlux,aeonFluxHits,aeonFluxQuality);
1265  miniMcPair->setDominatrack(aeonFlux);
1266  miniMcPair->setDominCommonHit(aeonFluxHits);
1267  miniMcPair->setAvgQuality(aeonFluxQuality);
1268  // special case for contam pairs
1269  //
1270  StContamPair* contamPair = dynamic_cast<StContamPair*>(miniMcPair);
1271  if(contamPair){
1272  const StMcTrack* mcTParent = (mcTrack ? mcTrack->parent() : 0);
1273 
1274  if (mcTParent){
1275  contamPair->setParentGeantId(mcTrack->parent()->geantId());
1276  contamPair->setPtMcParent(mcTrack->parent()->momentum().perp());
1277  contamPair->setEtaMcParent(mcTrack->parent()->pseudoRapidity());
1278  contamPair->setGeantProcess(mcTrack->startVertex()->geantProcess());
1279 
1280  contamPair->setStartX(mcTrack->startVertex()->position().x());
1281  contamPair->setStartY(mcTrack->startVertex()->position().y());
1282  contamPair->setStartZ(mcTrack->startVertex()->position().z());
1283 
1284  Int_t parentParentGeantId=0;
1285  Float_t parentParentPt=0;
1286  // check for parent of parent
1287  if(mcTrack->parent()->parent() &&
1288  mcTrack->parent()->parent()->geantId()>0){
1289  parentParentGeantId = mcTrack->parent()->parent()->geantId();
1290  parentParentPt = mcTrack->parent()->parent()->momentum().perp();
1291  }
1292  contamPair->setParentParentGeantId(parentParentGeantId);
1293  contamPair->setPtMcParentParent(parentParentPt);
1294 
1295 
1296  //* check if the parent doesnt start from the primary vertex
1297  /*
1298  if(mcTrack->parent()->parent() &&
1299  mcTrack->parent()->parent()->geantId()>0){
1300  cout << ">>WARNING: parent doesnt come from the primary vertex!" <<endl
1301  << "geant id : " << mcTrack->geantId() << endl
1302  << "r parent : " << mcTrack->parent()->startVertex()->position().perp() << endl
1303  << "parent geantId: " << mcTrack->parent()->geantId() << endl
1304  << "parent's parent: "
1305  << mcTrack->parent()->parent()->geantId() << endl
1306  << "pr rc pt: " << prTrack->geometry()->momentum().perp() << endl
1307  << "pr rc eta: " << prTrack->geometry()->momentum().pseudoRapidity() << endl
1308  << "rc fit pts: " << glTrack->fitTraits().numberOfFitPoints(kTpcId) << endl
1309  << "daughter pt: " << mcTrack->momentum().perp() << endl
1310  << "daughter p: " << mcTrack->momentum().mag() << endl
1311  << "daughter eta: " << mcTrack->pseudoRapidity() << endl
1312  << "parent pt: " << mcTrack->parent()->momentum().perp() << endl
1313  << "parent p: " << mcTrack->parent()->momentum().mag() << endl
1314  << "parent eta: " << mcTrack->parent()->pseudoRapidity() << endl
1315  << "parent's parent: " << mcTrack->parent()->parent()->momentum().perp() << endl;
1316  }
1317  */
1318  } else {
1319  StMiniMcMakerErrorCount++;
1320  if ( StMiniMcMakerErrorCount < 50){
1321  cout << "StMiniMcMaker::fillTrackPairInfo: WARNING mcTrack->parent() is NULL !! "
1322  << " If this comes from a normal simulation or embeding, please report !!"
1323  << " This is known to happen in pilup mode ONLY ... " << endl;
1324  }
1325  }
1326  }
1327 
1328 }
1329 
1330 //______________________________________________________________________________
1331 void StMiniMcMaker::fillRcTrackInfo(StTinyRcTrack* tinyRcTrack,
1332  const StTrack* prTrack,
1333  const StTrack* glTrack,
1334  Int_t nAssocMc)
1335 {
1336  if (!glTrack) {
1337  cout << "Error StMiniMcMaker::fillRcTrackInfo, glTrack pointer is zero " << glTrack << endl;
1338  return;
1339  }
1340  const StGlobalTrack *global = (StGlobalTrack*)glTrack;
1341  const StDcaGeometry* dcaGeo = global->dcaGeometry();
1342  tinyRcTrack->setValidGl();
1343  tinyRcTrack->setRecoKey(glTrack->key());
1344  tinyRcTrack->setDca00(1000.);
1345  if (dcaGeo ) {
1346  tinyRcTrack->setDca00(dcaGeo->impact());
1347  tinyRcTrack->setDca(1);
1348  StThreeVectorF glMom = dcaGeo->momentum();
1349  THelixTrack glHelix = dcaGeo->thelix();
1350 
1351 
1352  enum { kImpImp=0,
1353  kZImp, kZZ,
1354  kPsiImp, kPsiZ, kPsiPsi,
1355  kPtiImp, kPtiZ, kPtiPsi, kPtiPti,
1356  kTanImp, kTanZ, kTanPsi, kTanPti, kTanTan};
1357  const float *dcaErr=dcaGeo->errMatrix();
1358  // the indices of the error matrix correspond to
1359  // 0 - error on y (track position along pad row direction)
1360  // 1 - error on z (track position along drift direction)
1361  // 2 - error on Psi
1362  // 3 - error on pt,
1363  // 4 - error on tan(dipAngle)
1364 
1365  double pt = glMom.perp();
1366  tinyRcTrack->setPtGl(pt);
1367  tinyRcTrack->setPzGl(glMom.z());
1368  tinyRcTrack->setEtaGl(glMom.pseudoRapidity());
1369  tinyRcTrack->setPhiGl(glMom.phi());
1370  tinyRcTrack->setCurvGl(dcaGeo->curvature());
1371  tinyRcTrack->setTanLGl(dcaGeo->tanDip());
1372  tinyRcTrack->setSeedQuality(glTrack->seedQuality());
1373  float errorGl[5] =
1374  {float(dcaErr[kImpImp]) ,
1375  float(dcaErr[kZZ]) ,
1376  float(dcaErr[kPsiPsi]) ,
1377  float(dcaErr[kPtiPti]*pow(pt,4)) ,
1378  float(dcaErr[kTanTan])
1379  };
1380  for (int j=0;j<5;j++) {errorGl[j] = sqrt(errorGl[j]);}
1381  tinyRcTrack->setErrGl(errorGl);
1382  double vtx[3]={mRcVertexPos[0][0],mRcVertexPos[0][1],mRcVertexPos[0][2]};
1383  double dcaXY,dcaZ;
1384  glHelix.Dca(vtx,dcaXY,dcaZ,0);
1385  tinyRcTrack->setDcaXYGl(dcaXY);
1386  tinyRcTrack->setDcaZGl(dcaZ);
1387  tinyRcTrack->setDcaGl(sqrt(dcaXY*dcaXY+dcaZ*dcaZ));
1388  double mcv[3]={mMcVertexPos[0][0],mMcVertexPos[0][1],mMcVertexPos[0][2]};
1389  glHelix.Dca(mcv,dcaXY,dcaZ,0);
1390  tinyRcTrack->setDcaXYGlMcV(dcaXY);
1391  tinyRcTrack->setDcaZGlMcV(dcaZ);
1392  } else {
1393  tinyRcTrack->setDca(0);
1394  const StThreeVectorF& glMom = glTrack->geometry()->momentum();
1395 
1396  const StPhysicalHelixD& glHelix = glTrack->geometry()->helix();
1397 
1398 
1399 
1400  double pt = glMom.perp();
1401  tinyRcTrack->setPtGl(pt);
1402  tinyRcTrack->setPzGl(glMom.z());
1403  tinyRcTrack->setEtaGl(glMom.pseudoRapidity());
1404  tinyRcTrack->setPhiGl(glMom.phi());
1405  tinyRcTrack->setCurvGl(glTrack->geometry()->curvature());
1406  tinyRcTrack->setTanLGl(tan(glTrack->geometry()->dipAngle()));
1407  tinyRcTrack->setSeedQuality(glTrack->seedQuality());
1408  StMatrixF gCM = glTrack->fitTraits().covariantMatrix();
1409 //VP Float_t errorGl[5] = {gCM(1,1),gCM(2,2),gCM(3,3),gCM(4,4),gCM(5,5)};
1410  Float_t errorGl[5] = {
1411  Float_t(gCM[0][0]*pow(M_PI/180,2)) , //YY
1412  Float_t(gCM[1][1]) , //ZZ
1413  Float_t(gCM[3][3]*pow(M_PI/180,2)) , //PsiPsi
1414  Float_t(gCM[4][4]*pow(pt,4)) , //PtPt
1415  Float_t(gCM[2][2]) //tanLtanL
1416  };
1417  for (int j=0;j<5;j++) {errorGl[j] = sqrt(errorGl[j]);}
1418  tinyRcTrack->setErrGl(errorGl);
1419  //
1420  // reality check
1421  //
1422  // if(fabs(tinyRcTrack->mDcaPr - prTrack->impactParameter())>.001){
1423  // cout << " helix : " << tinyRcTrack->mDcaPr
1424  // << " impact: " << prTrack->impactParameter() << endl;
1425  // }
1426 
1427  tinyRcTrack->setDcaGl(glTrack->impactParameter());
1428  tinyRcTrack->setDcaXYGl(computeXY(mRcVertexPos,glTrack));
1429  //tinyRcTrack->setDcaZGl(computeZDca(mRcVertexPos,glTrack));
1430  tinyRcTrack->setDcaXYGl(computeXY(mRcVertexPos,glTrack));
1431  tinyRcTrack->setDcaZGlMcV(dcaz(glHelix,*mMcVertexPos,glTrack));
1432  tinyRcTrack->setDcaXYGlMcV(computeXY(mMcVertexPos,glTrack));
1433  }
1434 
1435  StDedxPidTraits* pid = findDedxPidTraits(glTrack);
1436  float meanDedx = (pid) ? pid->mean() : -999;
1437  tinyRcTrack->setDedx(meanDedx);
1438  short nDedxPts = (pid) ? pid->numberOfPoints() : 0;
1439  tinyRcTrack->setDedxPts(nDedxPts);
1440 
1441  //
1442  // common rc info
1443  //
1444  // first and last hit
1445  PAIRHIT hits = findFirstLastHit(glTrack);
1446  PAIRHIT fitHits = findFirstLastFitHit(glTrack);
1447 
1448  if(hits.first){
1449  // cout << "first hit z " << hits.first->position().z() << endl;
1450  tinyRcTrack->setFirstZ(hits.first->position().z());
1451  tinyRcTrack->setLastZ(hits.second->position().z());
1452  tinyRcTrack->setFirstPadrow(hits.first->padrow());
1453  tinyRcTrack->setLastPadrow(hits.second->padrow());
1454  tinyRcTrack->setFirstSector(hits.first->sector());
1455  tinyRcTrack->setLastSector(hits.second->sector());
1456  }
1457  else if ( glTrack->detectorInfo()->numberOfPoints()==0) {
1458  // could be an FTPC track, for which firstLastHit won't work
1459  cout << "Error: no hits?" << endl;
1460  cout << "Tpc points : " << glTrack->detectorInfo()->numberOfPoints(kTpcId) << endl;
1461  cout << "Svt points : " << glTrack->detectorInfo()->numberOfPoints(kTpcId) << endl;
1462  cout << "Ftpc points E: " << glTrack->detectorInfo()->numberOfPoints(kFtpcEastId) << endl;
1463  cout << "Ftpc points W: " << glTrack->detectorInfo()->numberOfPoints(kFtpcWestId) << endl;
1464  }
1465  if (fitHits.first) {
1466  tinyRcTrack->setFirstFitPadrow(fitHits.first->padrow());
1467  tinyRcTrack->setLastFitPadrow(fitHits.second->padrow());
1468  }
1469  else if ( glTrack->fitTraits().numberOfFitPoints()==0) {
1470  cout << "Error: no hit with usedInFit()>0" << endl;
1471  cout << "Tpc fit pts :" << glTrack->fitTraits().numberOfFitPoints(kTpcId) << endl;
1472  cout << "Svt fit pts :" << glTrack->fitTraits().numberOfFitPoints(kSvtId) << endl;
1473  cout << "Ftpc fit pts E:" << glTrack->fitTraits().numberOfFitPoints(kFtpcEastId) << endl;
1474  cout << "Ftpc fit pts C:" << glTrack->fitTraits().numberOfFitPoints(kFtpcWestId) << endl;
1475  }
1476 
1477  tinyRcTrack->setFitPts(glTrack->fitTraits().numberOfFitPoints(kTpcId));
1478  tinyRcTrack->setFitSvt(glTrack->fitTraits().numberOfFitPoints(kSvtId));
1479  tinyRcTrack->setFitSsd(glTrack->fitTraits().numberOfFitPoints(kSsdId));
1480  size_t ftpcFitPts = 0;
1481  if (tinyRcTrack->etaGl()>1.8)
1482  ftpcFitPts = glTrack->fitTraits().numberOfFitPoints(kFtpcWestId);
1483  if (tinyRcTrack->etaGl()<-1.8)
1484  ftpcFitPts = glTrack->fitTraits().numberOfFitPoints(kFtpcEastId);
1485  tinyRcTrack->setFitFtpc(ftpcFitPts);
1486  tinyRcTrack->setAllPts(glTrack->detectorInfo()->numberOfPoints(kTpcId));
1487  tinyRcTrack->setCharge(glTrack->geometry()->charge());
1488 
1489  tinyRcTrack->setNAssocMc(nAssocMc);
1490  tinyRcTrack->setNPossible(glTrack->numberOfPossiblePoints(kTpcId));
1491 
1492  if (prTrack) {
1493  tinyRcTrack->setValidPr();
1494  // with the introduction of the global track branch,
1495  // having the primary track pointer here is optional.
1496  const StThreeVectorF& prMom = prTrack->geometry()->momentum();
1497  const StPhysicalHelixD& prHelix = prTrack->geometry()->helix();
1498 
1499  double pt = prMom.perp();
1500  tinyRcTrack->setPtPr(pt);
1501  tinyRcTrack->setPzPr(prMom.z());
1502  tinyRcTrack->setEtaPr(prMom.pseudoRapidity());
1503  tinyRcTrack->setPhiPr(prMom.phi());
1504  tinyRcTrack->setCurvPr(prTrack->geometry()->curvature());
1505  tinyRcTrack->setTanLPr(tan(prTrack->geometry()->dipAngle()));
1506  tinyRcTrack->setChi2Pr(prTrack->fitTraits().chi2());
1507  tinyRcTrack->setFlag(prTrack->flag());
1508  tinyRcTrack->setDcaPr(prTrack->impactParameter());
1509  tinyRcTrack->setDcaXYPr(computeXY(mRcVertexPos,prTrack));
1510  //tinyRcTrack->setDcaZPr(computeZDca(mRcVertexPos,prTrack));
1511  tinyRcTrack->setDcaZPr(dcaz(prHelix,*mRcVertexPos));
1512  tinyRcTrack->setDcaXYPrMcV(computeXY(mMcVertexPos,prTrack));
1513  tinyRcTrack->setDcaZPrMcV(dcaz(prHelix,*mMcVertexPos));
1514  tinyRcTrack->setFitPts(prTrack->fitTraits().numberOfFitPoints(kTpcId));
1515  tinyRcTrack->setFitSvt(prTrack->fitTraits().numberOfFitPoints(kSvtId));
1516  tinyRcTrack->setFitSsd(prTrack->fitTraits().numberOfFitPoints(kSsdId));
1517  StMatrixF pCM = prTrack->fitTraits().covariantMatrix();
1518  Float_t errorPr[5] = {
1519  Float_t(pCM[0][0]*pow(M_PI/180,2)) , //YY
1520  Float_t(pCM[1][1]) , //ZZ
1521  Float_t(pCM[3][3]*pow(M_PI/180,2)) , //PsiPsi
1522  Float_t(pCM[4][4]*pow(pt,4)) , //PtPt
1523  Float_t(pCM[2][2]) //tanLtanL
1524  };
1525  for (int j=0;j<5;j++) {errorPr[j] = sqrt(errorPr[j]);}
1526  tinyRcTrack->setErrPr(errorPr);
1527  size_t ftpcFitPts = 0;
1528  if (tinyRcTrack->etaGl()>1.8)
1529  ftpcFitPts = prTrack->fitTraits().numberOfFitPoints(kFtpcWestId);
1530  if (tinyRcTrack->etaGl()<-1.8)
1531  ftpcFitPts = prTrack->fitTraits().numberOfFitPoints(kFtpcEastId);
1532  tinyRcTrack->setFitFtpc(ftpcFitPts);
1533 
1534  }
1535 
1536  //
1537  // EMC Information additions (Sept 2007) MCBS
1538  // Add to the track information about EMC towers.
1539  // - Find the tower at which the track points.
1540  // - Obtain all 8 towers (StEmcRawHit) surrounding the above tower (9 total, modulo edge effects)
1541  // - Sort all 9 towers according to their energy
1542  // - Store the ADC and energy of the 3 highest towers out of these 9
1543  // - Store the Soft Id of the highest tower
1544  // - Store the energy of StEmcPoint related to this track, if any.
1545  StEmcPosition emcPos;
1546  StThreeVectorD pos(0,0,0);
1547  StThreeVectorD mom(0,0,0);
1548  double magField = mRcEvent->runInfo()->magneticField();
1549 
1550  if (Debug()>2) {
1551  cout << "fillRcTrack, EMC information" << endl;
1552  cout << "Extrapolating to BEMC using B Field = " << magField*kilogauss/tesla << " tesla" << endl;
1553 // if (glTrack) {
1554 // cout << "hx curvature original from track: " << glTrack->outerGeometry()->curvature() << endl;
1555 // StPhysicalHelixD helix(glTrack->outerGeometry()->momentum(),glTrack->outerGeometry()->origin(),magField*kilogauss,glTrack->outerGeometry()->charge());
1556 // cout << "hx curvature as in EmcProjection: " << helix.curvature() << endl;
1557 // }
1558  }
1559 
1560  // Project Track onto BEMC radius.
1561  // Use the mag field obtained from run Info above, it
1562  // will get multiplied by tesla in StEmcPosition::projTrack.
1563  bool projOk;
1564 
1565  if (prTrack) {
1566  projOk = emcPos.trackOnEmc(&pos,&mom,prTrack,magField*kilogauss/tesla);
1567  }
1568  else {
1569  projOk = emcPos.trackOnEmc(&pos,&mom,glTrack,magField*kilogauss/tesla);
1570  }
1571  if (projOk) {
1572  // Track hits BEMC. Find the 9 closest towers around the projection.
1573  // Sort them according to energy. Keep the highest 3.
1574  int softIdProj(-1);
1575  std::vector<StEmcRawHit*> towersOfTrack;
1576  StEmcGeom* emcGeom = StEmcGeom::getEmcGeom("bemc");
1577  emcGeom->getId(pos.phi(),pos.pseudoRapidity(),softIdProj); // softId range [1,4800];
1578  for(int idEta=-1; idEta<2; ++idEta) {
1579  for (int idPhi=-1; idPhi<2; ++idPhi) {
1580  int towerId = emcPos.getNextTowerId(softIdProj,idEta,idPhi);
1581  if ( towerId == 0 )
1582  {
1583  continue; // off the end of the barrel
1584  }
1585  if (!emcGeom->checkId(towerId)) {
1586  // Valid tower SoftId. (Again, note that checkId returns "false" when the
1587  // softId is valid!)
1588  // Get StEmcRawHit corresponding to this softId
1589  // here is where the helper Emc container comes in handy
1590  StEmcRawHit* emcHit = mEmcIndex[towerId];
1591  if (emcHit!=0) {
1592  towersOfTrack.push_back(emcHit);
1593  } // check valid StEmcRawHit pointer
1594  }// check valid tower id
1595  }// dPhi loop
1596  }// dEta loop
1597  if (Debug()>1) {
1598  cout << "Outer Helix " << glTrack->outerGeometry()->helix() << endl;
1599  cout << "Track Projects to tower " << softIdProj << endl;
1600  cout << "Track hits at R= " << pos.perp() << " eta,phi: " << pos.pseudoRapidity() << ", " << pos.phi() << endl;
1601  cout << "Track Has " << towersOfTrack.size() << " total candidate towers" << endl;
1602  }
1603  if (towersOfTrack.size()>0) {
1604  // Obtained all towers (9 at most).
1605  // Sort them by energy, using helper function defined atop this file
1606  sort(towersOfTrack.begin(),towersOfTrack.end(),StEmcRawHitCompEne);
1607 
1608  // Store the ADC, Energy and SoftId of the 3 highest towers
1609  size_t maxTowers=3;
1610  if (towersOfTrack.size()<maxTowers) maxTowers=towersOfTrack.size();
1611  for (size_t iTow=0; iTow<maxTowers; ++iTow) {
1612  tinyRcTrack->setEmcTowerAdc(towersOfTrack[iTow]->adc(),iTow);
1613  tinyRcTrack->setEmcEnergyRcHit(towersOfTrack[iTow]->energy(),iTow);
1614  tinyRcTrack->setEmcSoftIdHiTowerRc(towersOfTrack[iTow]->softId(1),iTow);
1615  } //3 highest towers loop
1616 
1617  } // make sure there are towers
1618  else {
1619  for (size_t iTow=0; iTow<0; ++iTow) {
1620  tinyRcTrack->setEmcTowerAdc(-9,iTow);
1621  tinyRcTrack->setEmcEnergyRcHit(-9,iTow);
1622  tinyRcTrack->setEmcSoftIdHiTowerRc(-9,iTow);
1623  } //3 highest towers loop
1624 
1625  }
1626  if (Debug()>1) {
1627  cout << "rc track, key " << tinyRcTrack->recoKey() << endl;
1628  cout << "n Tpc Fit Pts " << tinyRcTrack->fitPts() << endl;
1629  cout << "3 bemc hit ene " << tinyRcTrack->emcEnergyRcHit(0)
1630  << ", " << tinyRcTrack->emcEnergyRcHit(1)
1631  << ", " << tinyRcTrack->emcEnergyRcHit(2) << endl;
1632  cout << "3 HiTow SoftId Rc " << tinyRcTrack->emcSoftIdHiTowerRc(0)
1633  << ", " << tinyRcTrack->emcSoftIdHiTowerRc(1)
1634  << ", " << tinyRcTrack->emcSoftIdHiTowerRc(2) << endl;
1635  float etaTow, phiTow;
1636  emcGeom->getEtaPhi(tinyRcTrack->emcSoftIdHiTowerRc(0),etaTow,phiTow);
1637  cout << "Hi Tow eta, phi Rc " << etaTow << ", " << phiTow << endl;
1638 
1639  }
1640  }// track has valid projection to EMC.
1641  return;
1642 }
1643 
1644 //______________________________________________________________________________
1645 void StMiniMcMaker::fillMcTrackInfo(StTinyMcTrack* tinyMcTrack,
1646  const StMcTrack* mcTrack,
1647  Int_t nAssocGl, Int_t nAssocPr)
1648 {
1649  if (mcTrack) {
1650  tinyMcTrack->setValid();
1651  const StThreeVectorF& mcMom = mcTrack->momentum();
1652 
1653  tinyMcTrack->setKey(mcTrack->key());
1654  tinyMcTrack->setPrimary(mcTrack->IsPrimary());
1655  tinyMcTrack->setPtMc(mcMom.perp());
1656  tinyMcTrack->setPzMc(mcMom.z());
1657  tinyMcTrack->setEtaMc(mcMom.pseudoRapidity());
1658  tinyMcTrack->setPhiMc(mcMom.phi());
1659  tinyMcTrack->setNHitMc(mcTrack->tpcHits().size());
1660  tinyMcTrack->setNSvtHitMc(mcTrack->svtHits().size());
1661  tinyMcTrack->setNSsdHitMc(mcTrack->ssdHits().size());
1662  tinyMcTrack->setNFtpcHitMc(mcTrack->ftpcHits().size());
1663  tinyMcTrack->setGeantId(mcTrack->geantId());
1664  tinyMcTrack->setPdgId(mcTrack->pdgId());
1665  short chargeMc = -9999;
1666  if (mcTrack->particleDefinition()) chargeMc = static_cast<short>(mcTrack->particleDefinition()->charge());
1667  tinyMcTrack->setChargeMc(chargeMc);
1668 
1669  tinyMcTrack->setNAssocGl(nAssocGl);
1670  tinyMcTrack->setNAssocPr(nAssocPr);
1671 
1672  float stopR=(mcTrack->stopVertex()) ? mcTrack->stopVertex()->position().perp() : 999;
1673  tinyMcTrack->setStopR(stopR);
1674 
1675  // if(stopR<999) cout << ">>stop r=" << stopR << endl;
1676  tinyMcTrack->setKey(mcTrack->key());
1677  if (mcTrack->parent()!=0) {
1678  tinyMcTrack->setParentKey(mcTrack->parent()->key());
1679  tinyMcTrack->setParentGeantId(mcTrack->parent()->geantId());
1680  }
1681  tinyMcTrack->setPtMc(mcMom.perp());
1682  tinyMcTrack->setPzMc(mcMom.z());
1683  tinyMcTrack->setEtaMc(mcMom.pseudoRapidity());
1684  tinyMcTrack->setPhiMc(mcMom.phi());
1685  tinyMcTrack->setNHitMc(mcTrack->tpcHits().size());
1686  tinyMcTrack->setNSvtHitMc(mcTrack->svtHits().size());
1687  tinyMcTrack->setNFtpcHitMc(mcTrack->ftpcHits().size());
1688  tinyMcTrack->setNBemcHitMc(mcTrack->bemcHits().size());
1689  tinyMcTrack->setNBprsHitMc(mcTrack->bprsHits().size());
1690  tinyMcTrack->setNBsmdeHitMc(mcTrack->bsmdeHits().size());
1691  tinyMcTrack->setNBsmdpHitMc(mcTrack->bsmdpHits().size());
1692  tinyMcTrack->setNEemcHitMc(mcTrack->eemcHits().size());
1693  tinyMcTrack->setNBprsHitMc(mcTrack->bprsHits().size());
1694  tinyMcTrack->setNBsmdeHitMc(mcTrack->bsmdeHits().size());
1695  tinyMcTrack->setNBsmdpHitMc(mcTrack->bsmdpHits().size());
1696  tinyMcTrack->setGeantId(mcTrack->geantId());
1697 
1698  tinyMcTrack->setNAssocGl(nAssocGl);
1699  tinyMcTrack->setNAssocPr(nAssocPr);
1700 
1701  // Fill Mc Emc constMcTrack Information
1702  // Find the 3 highest energy towers and store them
1703  // This is done by copying the hits to a vector for sorting
1704  // according to the hit energy (dE). Then the first 3 values
1705  // i.e. the highest values, of dE are copied into the StTinyMcTrack
1706  // The total sum of dE for all Calorimeter hits is also computed and stored.
1707  // The softId of the hit with the highest tower is also stored.
1708  StPtrVecMcCalorimeterHit bemcHits = mcTrack->bemcHits();
1709  if (bemcHits.size()>0) {
1710  std::vector<StMcCalorimeterHit*> bemcHitsSorted(bemcHits.size());
1711  int hiEnergyHitSoftId(-999);
1712  double sumEnergy(0);
1713  copy(bemcHits.begin(),bemcHits.end(),bemcHitsSorted.begin());
1714  sort(bemcHitsSorted.begin(),bemcHitsSorted.end(),StMcCalorimeterHitCompdE);
1715  if (Debug()>2) {
1716  cout << "fillMcTrackInfo() Check sorting of dE for StMcCalorimeterHits" << endl;
1717  }
1718  StEmcGeom* emcGeom = StEmcGeom::getEmcGeom("bemc");
1719  for (StMcCalorimeterHitIterator bhi=bemcHitsSorted.begin();
1720  bhi!=bemcHitsSorted.end();
1721  ++bhi) {
1722  StMcCalorimeterHit* bh = *bhi;
1723  float eta(0);
1724  emcGeom->getEta(bh->module(),bh->eta(),eta);
1725  sumEnergy += (*bhi)->dE()*scaleFactor(eta);
1726  if (Debug()>2) cout << bh->dE()*scaleFactor(eta) << endl;
1727  }
1728  // Fill top 3 towers into track class (energy and SoftId), but only if there are enough hits.
1729  size_t maxHits = 3;
1730  if (bemcHitsSorted.size()<maxHits) maxHits=bemcHitsSorted.size();
1731  for (size_t iCalHit=0; iCalHit<maxHits; ++iCalHit) {
1732  float eta(0);
1733  emcGeom->getEta(bemcHitsSorted[iCalHit]->module(),bemcHitsSorted[iCalHit]->eta(),eta);
1734  tinyMcTrack->setEmcEnergyMcHit(bemcHitsSorted[iCalHit]->dE()*scaleFactor(eta),iCalHit);
1735  emcGeom->getId(bemcHitsSorted[iCalHit]->module(),
1736  bemcHitsSorted[iCalHit]->eta(),
1737  bemcHitsSorted[iCalHit]->sub(),hiEnergyHitSoftId);
1738  tinyMcTrack->setEmcSoftIdHiTowerMc(static_cast<Short_t>(hiEnergyHitSoftId),iCalHit);
1739  }
1740  // Fill the rest of the hits, from maxHits up to 3, with zeros
1741  for (size_t iCalHit2=maxHits; iCalHit2<3; ++iCalHit2) {
1742  tinyMcTrack->setEmcEnergyMcHit(-9,iCalHit2);
1743  tinyMcTrack->setEmcSoftIdHiTowerMc(-9,iCalHit2);
1744 
1745  }
1746  tinyMcTrack->setEmcEnergyMcSum(sumEnergy);
1747 
1748  if (Debug()>1) {
1749  cout << "mc track, geant Id " << tinyMcTrack->geantId() << endl;
1750  cout << "parent geantId " << tinyMcTrack->parentGeantId() << endl;
1751  cout << "n Tpc Hits " << tinyMcTrack->nHitMc() << endl;
1752  cout << "n Bemc Hits " << tinyMcTrack->nBemcHitMc() << endl;
1753  cout << "3 bemc hit ene " << tinyMcTrack->emcEnergyMcHit(0)
1754  << ", " << tinyMcTrack->emcEnergyMcHit(1)
1755  << ", " << tinyMcTrack->emcEnergyMcHit(2) << endl;
1756  cout << "MC 3 Hi SoftIds " << tinyMcTrack->emcSoftIdHiTowerMc(0)
1757  << ", " << tinyMcTrack->emcSoftIdHiTowerMc(1)
1758  << ", " << tinyMcTrack->emcSoftIdHiTowerMc(2) << endl;
1759  cout << "emc energy sum " << tinyMcTrack->emcEnergyMcSum() << endl;
1760  cout << "MC trk momentum " << tinyMcTrack->pMc() << endl;
1761  cout << "MC eta, phi " << tinyMcTrack->etaMc() << ", " << tinyMcTrack->phiMc() << endl;
1762  float etaTow, phiTow;
1763  emcGeom->getEtaPhi(tinyMcTrack->emcSoftIdHiTowerMc(0),etaTow,phiTow);
1764  cout << "MC HiTow eta,phi " << etaTow << ", " << phiTow << endl;
1765  if (mEmcIndex[tinyMcTrack->emcSoftIdHiTowerMc(0)]) {
1766  StEmcRawHit* rawHit = mEmcIndex[tinyMcTrack->emcSoftIdHiTowerMc(0)];
1767  cout << "RC Ene for id " << rawHit->energy()
1768  << " m= " << rawHit->module()
1769  << " e= " << rawHit->eta()
1770  << " s= " << rawHit->sub()
1771  << endl;
1772  //cout << "RC id rawHit " << rawHit->softId(kBarrelEmcTowerId) << endl; //doesn't work?
1773  }
1774  else {
1775  cout << "Soft Id Mc " << tinyMcTrack->emcSoftIdHiTowerMc(0) << " has no StEmcRawHit" << endl;
1776  }
1777  } // debug
1778 
1779  } // if the track has bemc hits
1780  // if(stopR<999) cout << ">>stop r=" << stopR << endl;
1781  }// if (mcTrack)
1782  return;
1783 }
1784 /*
1785  given a mc track, returns a vector of matched associated pairs.
1786 */
1787 //______________________________________________________________________________
1788 StTrackPairInfo* StMiniMcMaker::findBestMatchedGlobal(const StMcTrack* mcTrack)
1789 {
1790  pair<mcTrackMapIter,mcTrackMapIter> mcBounds = mMcTrackMap->equal_range((const StMcTrack*)mcTrack);
1791  StTrackPairInfo* candTrackPair = 0; // used for finding the best matched track
1792  const StGlobalTrack* candTrack = 0;
1793  mcTrackMapIter mcMapIter = mcBounds.first;
1794  for ( ; mcMapIter != mcBounds.second; ++mcMapIter){
1795  StTrackPairInfo* assocPair = (*mcMapIter).second;
1796  const StGlobalTrack* globTrack = assocPair->partnerTrack();
1797  if (Debug() > 1) {
1798  cout << * assocPair << endl;
1799  cout << "globTrack FitPoints Tpc/FtpcE/W = " << globTrack->fitTraits().numberOfFitPoints(kTpcId)
1800  << "/" << globTrack->fitTraits().numberOfFitPoints(kFtpcEastId)
1801  << "/" << globTrack->fitTraits().numberOfFitPoints(kFtpcWestId) << endl;
1802  }
1803  if (!globTrack || globTrack->flag()<=0) continue;
1804  if (globTrack->fitTraits().numberOfFitPoints(kTpcId)>=10 ||
1805  globTrack->fitTraits().numberOfFitPoints(kFtpcEastId)>=5 ||
1806  globTrack->fitTraits().numberOfFitPoints(kFtpcWestId)>=5) {
1807  if (!candTrackPair) {
1808  candTrackPair = assocPair;
1809  candTrack = globTrack;
1810  }
1811  else if (globTrack->fitTraits().numberOfFitPoints(kTpcId) > candTrack->fitTraits().numberOfFitPoints(kTpcId)) {
1812  candTrackPair = assocPair;
1813  candTrack = globTrack;
1814  }
1815  else if (globTrack->fitTraits().numberOfFitPoints(kFtpcEastId) > candTrack->fitTraits().numberOfFitPoints(kFtpcEastId)) {
1816  candTrackPair = assocPair;
1817  candTrack = globTrack;
1818  }
1819  else if (globTrack->fitTraits().numberOfFitPoints(kFtpcWestId) > candTrack->fitTraits().numberOfFitPoints(kFtpcWestId)) {
1820  candTrackPair = assocPair;
1821  candTrack = globTrack;
1822  }
1823 
1824  } // fit points requirement
1825  }// bounds loop
1826  return candTrackPair; // Note that candTrack might be zero, for example if only one track is matched and has 9 tpc fit pts.
1827 }
1828 //______________________________________________________________________________
1829 PAIRVEC StMiniMcMaker::findMatchedRc(const StMcTrack* mcTrack)
1830 {
1831  //
1832  // now find the associated tracks
1833  //
1834  pair<mcTrackMapIter,mcTrackMapIter> mcBounds = mMcTrackMap->equal_range((const StMcTrack*)mcTrack);
1835 
1836  PAIRVEC candPair; // used for finding the best matched track
1837 
1838  mcTrackMapIter mcMapIter = mcBounds.first;
1839  for( ; mcMapIter != mcBounds.second; mcMapIter++){
1840  StTrackPairInfo* assocPair = (*mcMapIter).second;
1841 
1842  const StGlobalTrack* glTrack = assocPair->partnerTrack();
1843 
1844  //
1845  // primary tracks only
1846  //
1847  const StPrimaryTrack* prTrack=0;
1848  if(!(prTrack = isPrimaryTrack(glTrack))) continue;
1849 
1850  // conspicuously absent is a pt cut
1851 
1852  //
1853  // share enough hits?
1854  //
1855  // if(assocPair->percentOfPairedTpcHits()<mSharedHitsCut) continue;
1856 
1857  //
1858  // bare minimum cut (fit hits, etc.)
1859  //
1860  if(!accept(glTrack) || !accept(prTrack)) continue;
1861 
1862  //
1863  // same sign?
1864  //
1865  //if(!isSameSign(prTrack,mcTrack)) continue;
1866 
1867  //
1868  // ok, everything's ok. save the global track (which is a primary trk).
1869  // we actually need both the global and primary track info
1870  //
1871  candPair.push_back(assocPair);
1872 
1873  } // map iter
1874 
1875  return candPair;
1876 }
1877 
1878 //______________________________________________________________________________
1879 PAIRHIT StMiniMcMaker::findFirstLastHit(const StTrack* track)
1880 {
1881  const StPtrVecHit& hits = track->detectorInfo()->hits(kTpcId);
1882  std::vector<const StTpcHit*> vec;
1883  // fill the vector
1884  // cout << "\thits size " << hits.size() << endl;
1885  for(UInt_t i=0; i<hits.size(); i++){
1886  const StTpcHit* hit = dynamic_cast<const StTpcHit*>(hits[i]);
1887  if(!hit) continue;
1888  vec.push_back(hit);
1889  }
1890  sort(vec.begin(),vec.end(),hitCmp);
1891  if (vec.size()) {
1892  return PAIRHIT(vec[0],vec[vec.size()-1]);
1893  }
1894  else {
1895  const StTpcHit* empty = 0;
1896  return PAIRHIT(empty,empty);
1897  }
1898 }
1899 
1900 //______________________________________________________________________________
1901 PAIRHIT StMiniMcMaker::findFirstLastFitHit(const StTrack* track)
1902 {
1903  const StPtrVecHit& hits = track->detectorInfo()->hits(kTpcId);
1904  std::vector<const StTpcHit*> vec;
1905  // fill the vector
1906 
1907  for(UInt_t i=0; i<hits.size(); i++){
1908  const StTpcHit* hit = dynamic_cast<const StTpcHit*>(hits[i]);
1909  if(!hit) continue;
1910  if(!hit->usedInFit()) continue;
1911  vec.push_back(hit);
1912  }
1913  sort(vec.begin(),vec.end(),hitCmp);
1914  if (vec.size()) {
1915  return PAIRHIT(vec[0],vec[vec.size()-1]);
1916  }
1917  else {
1918  const StTpcHit* empty = 0;
1919  return PAIRHIT(empty,empty);
1920  }
1921 }
1922 
1923 //--------- SIMPLE HELPERS--------------
1924 
1925 //______________________________________________________________________________
1926 // bare minimum cuts to accept a raw mc track
1927 Bool_t StMiniMcMaker::acceptRaw(const StMcTrack* mcTrack)
1928 {
1929  return (1);
1930  // commented out the cuts below, as they made sense for
1931  // single particle spectra in the TPC, but minimc is now
1932  // more general, so as long as there is a valid track pointer
1933  // we take it.
1934  //return (mcTrack &&
1935  // mcTrack->particleDefinition()->charge()!=0&&
1936  //fabs(mcTrack->momentum().pseudoRapidity())<=4.);
1937 }
1938 
1939 //______________________________________________________________________________
1940 // cut on an mc track to determine if we should look
1941 // for a matched rc track.
1942 Bool_t StMiniMcMaker::accept(const StMcTrack* mcTrack)
1943 {
1944  return (mcTrack && mcTrack->tpcHits().size() >= 10);
1945 }
1946 
1947 //______________________________________________________________________________
1948 // bare minimum cut to accept a rc track as valid
1949 Bool_t StMiniMcMaker::accept(const StTrack* rcTrack)
1950 {
1951  UInt_t nFitPoint = rcTrack->fitTraits().numberOfFitPoints(kUnknownId);
1952  return ( nFitPoint>=5 && rcTrack->flag()>0 );
1953 }
1954 
1955 //______________________________________________________________________________
1956 // are the mc and rc track of the same sign?
1957 Bool_t StMiniMcMaker::isSameSign(const StTrack* rcTrack,const StMcTrack* mcTrack)
1958 {
1959  return (rcTrack->geometry()->charge()*
1960  mcTrack->particleDefinition()->charge()>0);
1961 }
1962 
1963 //______________________________________________________________________________
1964 // possible pt cut when looking for split, background rc tracks
1965 Bool_t StMiniMcMaker::acceptPt(const StTrack* track)
1966 {
1967  return (track->geometry()->momentum().perp()>=mMinPt &&
1968  track->geometry()->momentum().perp()<=mMaxPt);
1969 }
1970 
1971 //______________________________________________________________________________
1972 // possible pt cut on mc tracks
1973 Bool_t StMiniMcMaker::acceptPt(const StMcTrack *track)
1974 {
1975  return (track->momentum().perp()>=mMinPt &&
1976  track->momentum().perp()<=mMaxPt);
1977 }
1978 //______________________________________________________________________________
1979 // cut for debugging purposes
1980 Bool_t StMiniMcMaker::acceptDebug(const StMcTrack* track)
1981 {
1982  return (track->momentum().pseudoRapidity() <= .1
1983  && track->momentum().pseudoRapidity() >= 0 &&
1984  (track->geantId()==9 || track->geantId()==12 || track->geantId()==15));
1985 
1986 }
1987 
1988 //______________________________________________________________________________
1989 // cut for flow centrality definition
1990 Bool_t StMiniMcMaker::acceptCentrality(const StTrack* track)
1991 {
1992  return (fabs(track->geometry()->momentum().pseudoRapidity())<.75);
1993 }
1994 
1995 //______________________________________________________________________________
1996 // cut for h- uncorrected centrality
1997 Bool_t StMiniMcMaker::acceptUncorrected(const StTrack* track)
1998 {
1999  return (
2000  track->geometry()->charge()<0 &&
2001  track->fitTraits().numberOfFitPoints(kTpcId)>=10 &&
2002  fabs(track->geometry()->momentum().pseudoRapidity())<=0.5 &&
2003  track->geometry()->helix().distance(*mRcVertexPos)<3
2004  );
2005 }
2006 
2007 //______________________________________________________________________________
2008 Bool_t StMiniMcMaker::acceptGlobals(const StTrack* track)
2009 {
2010  return 1;
2011 }
2012 
2013 //______________________________________________________________________________
2014 Bool_t StMiniMcMaker::acceptFTPC(const StTrack* prTrack)
2015 {
2016  return 1;
2017 }
2018 
2019 //______________________________________________________________________________
2020 // positive track flag
2021 Bool_t StMiniMcMaker::ok(const StTrack* track)
2022 {
2023  return (track && track->flag()>0);
2024 }
2025 //______________________________________________________________________________
2026 // Good Global RC Tracks with fitpts > 20
2027 Bool_t StMiniMcMaker::acceptGood20(const StTrack* track)
2028 {
2029  UInt_t nFitPoint = track->fitTraits().numberOfFitPoints(kTpcId);
2030  return (track && nFitPoint >= 20);
2031 }
2032 
2033 //______________________________________________________________________________
2034 // Good Global MC Tracks with fitpts > 20
2035 Bool_t StMiniMcMaker::acceptGood20(const StMcTrack* track)
2036 {
2037  return (track && track->tpcHits().size() >= 20);
2038 }
2039 
2040 //______________________________________________________________________________
2041 // checks if the rc track is from the vertex
2042 const StPrimaryTrack* StMiniMcMaker::isPrimaryTrack(const StTrack* glTrack)
2043 {
2044  if(!glTrack) return 0;
2045  return dynamic_cast<const StPrimaryTrack*>(glTrack->node()->track(primary));
2046 }
2047 
2048 //______________________________________________________________________________
2049 // checks if the mc track is from the vertex
2050 Bool_t StMiniMcMaker::isPrimaryTrack(const StMcTrack* mcTrack)
2051 {
2052 
2053  return(mcTrack->startVertex() == mMcEvent->primaryVertex());
2054 
2055 }
2056 
2057 //______________________________________________________________________________
2058 // xy dca
2059 Float_t StMiniMcMaker::computeXY(const StThreeVectorF* pos, const StTrack* track)
2060 {
2061  //
2062  // find the distance between the center of the circle and pos.
2063  // if the radius of curvature > distance, then call
2064  // it positive.
2065  //
2066  double xCenter = track->geometry()->helix().xcenter();
2067  double yCenter = track->geometry()->helix().ycenter();
2068  double radius = 1.0/track->geometry()->helix().curvature();
2069 
2070  double dPosCenter
2071  = TMath::Sqrt( (pos->x()-xCenter) * (pos->x()-xCenter) +
2072  (pos->y()-yCenter) * (pos->y()-yCenter));
2073 
2074  return (Float_t) ( radius - dPosCenter );
2075 }
2076 
2077 //______________________________________________________________________________
2078 // z dca from ben norman (no longer used)
2079 Float_t StMiniMcMaker::computeZDca(const StThreeVectorF* point, const StTrack* track)
2080 {
2081  const StPhysicalHelixD& helix = track->geometry()->helix();
2082  pairD path = helix.pathLength(point->perp());
2083 
2084  const StThreeVectorD& pos1 = helix.at(path.first);
2085  const StThreeVectorD& pos2 = helix.at(path.second);
2086  const StThreeVectorD dis1 = *point - pos1;
2087  const StThreeVectorD dis2 = *point - pos2;
2088 
2089  double dcaZ = (dis1.mag() < dis2.mag()) ? dis1.z() : dis2.z();
2090  if(std::isnan(dcaZ)) return 999;
2091  return dcaZ;
2092 }
2093 //______________________________________________________________________________
2094 StDedxPidTraits* StMiniMcMaker::findDedxPidTraits(const StTrack* track)
2095 {
2096  StDedxPidTraits* pid=0;
2097  StPtrVecTrackPidTraits traits = track->pidTraits(kTpcId);
2098 
2099  for (UInt_t i = 0; i < traits.size(); i++) {
2100  pid = dynamic_cast<StDedxPidTraits*>(traits[i]);
2101  if (pid && pid->method() == kLikelihoodFitId) break;
2102  }
2103  return pid;
2104 }
2105 
2106 //______________________________________________________________________________
2107 // debugging
2108 void StMiniMcMaker::checkMerged(const StMcTrack* mergedMcTrack, Int_t mergedCommonHits,
2109  const StTrack* track)
2110 {
2111 
2112  if(!isPrimaryTrack(mergedMcTrack))
2113  cout << "NOT a MC PRIMARY track" << endl;
2114 
2115  cout << "rc key : " << track->key() << endl;
2116  cout << "rc pr pt : "
2117  << track->geometry()->momentum().perp() << endl;
2118  cout << "rc pr flag : " << track->flag() << endl;
2119  cout << "rc pr dca : "
2120  << track->geometry()->helix().distance(*mRcVertexPos) << endl;
2121  cout << "mc key : " << mergedMcTrack->key() << endl;
2122  cout << "mc evt gen : " << mergedMcTrack->eventGenLabel() << endl;
2123  cout << "mc pt : "
2124  << mergedMcTrack->momentum().perp() << endl;
2125  cout << "mc eta: " << mergedMcTrack->momentum().pseudoRapidity() << endl;
2126  cout << "mc hits : "
2127  << mergedMcTrack->tpcHits().size() << endl;
2128  cout << "geant : " << mergedMcTrack->geantId() << endl;
2129  cout << "common hits : "
2130  << mergedCommonHits << endl;
2131  cout << "rc all hits : "
2132  << track->detectorInfo()->numberOfPoints(kTpcId) << endl;
2133  cout << "rc fit hits : "
2134  << track->fitTraits().numberOfFitPoints(kTpcId) << endl;
2135 
2136 }
2137 
2138 //______________________________________________________________________________
2139 // debugging
2140 void StMiniMcMaker::checkSplit(const StMcTrack* mcTrack, const StTrack* glTrack,
2141  Int_t commonHits)
2142 {
2143  // find the best rc match to this mc track
2144  //
2145  PAIRVEC testPair = findMatchedRc(mcTrack);
2146  sort(testPair.begin(), testPair.end(), sortCmp);
2147 
2148  cout << "#######################################" << endl;
2149  cout << "CHECK SPLIT:" << endl;
2150  cout << "mc pt : "
2151  << mcTrack->momentum().perp() << endl;
2152  cout << "mc key: "
2153  << mcTrack->key() << endl;
2154  if (mcTrack->particleDefinition()) {
2155  cout << "mc charge: "
2156  << mcTrack->particleDefinition()->charge() << endl;
2157  }
2158  cout << "mc eta : "
2159  << mcTrack->momentum().pseudoRapidity() << endl;
2160  cout << "mc pts: "
2161  << mcTrack->tpcHits().size() << endl;
2162  cout << "gl pt : "
2163  << glTrack->geometry()->momentum().perp() << endl;
2164  cout << "all hits : "
2165  << glTrack->detectorInfo()->numberOfPoints(kTpcId) << endl;
2166  cout << "fit pts : "
2167  << glTrack->fitTraits().numberOfFitPoints(kTpcId)<< endl;
2168  cout << "common hits : " << commonHits << endl;
2169 
2170  // if(mcTrack->startVertex() == mMcEvent->primaryVertex()){
2171  // cout << "YES mc track is a primary" << endl;
2172  // }
2173  if(testPair.size() == 1) {
2174  cout << "ERROR in split. not really a split track. "
2175  << "maybe due to pt cut?" << endl;
2176  cout << "mc start vertex : " << mcTrack->startVertex()->position() << endl;
2177 
2178  }
2179  cout << ">>Here are the rc tracks" << endl;
2180  for(unsigned int i=0; i<testPair.size(); i++){
2181  cout << "i: " << i << endl;
2182  const StGlobalTrack* testGlTrack = testPair[i]->partnerTrack();
2183  cout << "pt rc : "
2184  << testGlTrack->geometry()->momentum().perp() << endl;
2185  cout << "all hits : "
2186  << testGlTrack-> detectorInfo()->numberOfPoints(kTpcId) << endl;
2187  cout << "fit pts : "
2188  << testGlTrack->fitTraits().numberOfFitPoints(kTpcId)<< endl;
2189  cout << "common hits : "
2190  << testPair[i]->commonTpcHits() << endl;
2191  }
2192  cout << "#######################################" << endl;
2193 
2194 }
2195 
2196 
2197 //______________________________________________________________________________
2198 void StMiniMcMaker::checkContam(const StMcTrack* mcTrack, const StGlobalTrack* glTrack,
2199  Int_t commonHits)
2200 {
2201 
2202  cout << "############## " << endl;
2203  cout << "CHECK CONTAM " << endl;
2204 
2205  if(mcTrack->startVertex() == mMcEvent->primaryVertex()){
2206  cout << "\tERROR mc track is a primary!\n" << endl;
2207  }
2208  const StPrimaryTrack* prTrack=isPrimaryTrack(glTrack);
2209  cout << "common hits : " << commonHits << endl
2210  << "rc key : " << prTrack->key() << endl
2211  << "mc key : " << mcTrack->key() << endl
2212  << "pr pt : " << prTrack->geometry()->momentum().perp() << endl
2213  << "mc pt : " << mcTrack->momentum().perp() << endl
2214  << "fit pts : " << glTrack->fitTraits().numberOfFitPoints(kTpcId)
2215  << endl
2216  << "all pts : "<< glTrack->detectorInfo()->numberOfPoints(kTpcId)<< endl
2217  << "gl dca xy : " << computeXY(mRcVertexPos,glTrack) << endl;
2218  cout << ">>Here are all the mc tracks matched to this rc track" << endl;
2219 
2220  pair<rcTrackMapIter,rcTrackMapIter> rcBounds
2221  = mRcTrackMap->equal_range(glTrack);
2222 
2223  rcTrackMapIter rcMapIter = rcBounds.first;
2224 
2225  for( ; rcMapIter != rcBounds.second; rcMapIter++){
2226  StTrackPairInfo* assocPair = (*rcMapIter).second;
2227  const StMcTrack* mcCandTrack = assocPair->partnerMcTrack();
2228 
2229  cout << "common hits=" << assocPair->commonTpcHits()
2230  << ", is mc primary=" << (isPrimaryTrack(mcTrack)?"yes":"no")
2231  << ", mc key=" << mcCandTrack->key()
2232  << ", hits=" << mcCandTrack->tpcHits().size() << endl;
2233  }
2234  cout << ">>Here are the rc tracks matched to this mc track" << endl;
2235 
2236  PAIRVEC testPair = findMatchedRc(mcTrack);
2237  sort(testPair.begin(), testPair.end(), sortCmp);
2238 
2239  for(unsigned int i=0; i<testPair.size(); i++){
2240  cout << "i: " << i << endl;
2241  const StGlobalTrack* testGlTrack = testPair[i]->partnerTrack();
2242 
2243  cout << "common hits=" << testPair[i]->commonTpcHits()
2244  << ", rc key=" << testGlTrack->key() << endl;
2245  }
2246 
2247 }
2248 
2249 //______________________________________________________________________________
2250 size_t StMiniMcMaker::getIndex(size_t mult)
2251 {
2252 
2253  // note: this depends on the production
2254  //
2255 
2256  // P02gd, 2k2 data, Nch cuts
2257  if (mult >= 510) return 0;
2258  if (mult >= 431) return 1;
2259  if (mult >= 312) return 2;
2260  if (mult >= 217) return 3;
2261  if (mult >= 146) return 4;
2262  if (mult >= 94 ) return 5;
2263  if (mult >= 56 ) return 6;
2264  if (mult >= 30 ) return 7;
2265  if (mult >= 14 ) return 8;
2266  return 9;
2267 }
2268 
2269 //______________________________________________________________________________
2270 void StMiniMcMaker::AppendMCDaughterTrack()
2271 {
2272  if (Debug())
2273  cout << "##StMiniMcMaker::AppendMCDaughterTrack()"<< endl;
2274 
2275 // std::vector<int> enteredGlobalTracks;
2276  const StPtrVecMcTrack& allmcTracks = mMcEvent->tracks();
2277  cout << "size of mcEvent->tracks() : "<< allmcTracks.size() << endl;
2278 
2279  Int_t nAppendMC(0);
2280  StMcTrackConstIterator allMcTrkIter = allmcTracks.begin();
2281  for (; allMcTrkIter != allmcTracks.end(); ++allMcTrkIter) {
2282  const StMcTrack* mcGlobTrack = *allMcTrkIter;
2283  if (isPrimaryTrack(mcGlobTrack)) continue;
2284  if (!acceptRaw(mcGlobTrack)) continue; // loose eta cut (4 units, so should include ftpc).
2285 // if (accept(mcGlobTrack) || mcGlobTrack->ftpcHits().size()>=5) { // 10 tpc hits or 5 ftpc hits
2286  if (Debug()>1)
2287  cout << "accepted mc global track, key "<< mcGlobTrack->key() << endl;
2288  // Ok, track is accepted, query the map for its reco tracks.
2289 
2290  StTinyMcTrack tinyMcTrack;
2291  fillMcTrackInfo(&tinyMcTrack, mcGlobTrack, 0, 0);
2292  mMiniMcEvent->addMcTrack(&tinyMcTrack);
2293  nAppendMC++;
2294 
2295 // } // mc hits condition
2296  } // end of global track match loop
2297 
2298  cout << "\tappended mc tracks: "<< nAppendMC << endl;
2299 }
2300 
2301 //______________________________________________________________________________
2302 void StMiniMcMaker::dominatTkInfo(const StTrack* recTrack,int &dominatrackKey ,int& dominatrackHits,float& avgQuality) {
2303  // initialize return values.
2304  // dominatrack key initialized to nonsense, quality initialized to 0.
2305  // Note, I'm using shorts, which should be ok up to 32768, but if we
2306  // ever have track keys above this in an event, there will be trouble.
2307 
2308  int DetectorList[kMaxDetectorId]={0};
2309  typedef std::map< int,float> myMap_t;
2310  typedef myMap_t::const_iterator myIter_t;
2311  myMap_t idTruths;
2312 
2313  const StPtrVecHit &recHits = recTrack->detectorInfo()->hits();
2314 // Loop to store all the mc track keys and quality of every reco hit on the track.
2315  int nHits = recHits.size();
2316  for (int hi=0;hi<nHits; hi++) {
2317  const StHit* rHit = recHits[hi];
2318  int id = rHit->idTruth(); if (!id) continue;
2319  int qa = rHit->qaTruth(); if (!qa) qa = 1;
2320  idTruths[id]+=qa;
2321  }
2322  int tkBest=0; float qaBest=0,qaSum=0;
2323  for (myIter_t it=idTruths.begin(); it!=idTruths.end();++it) {
2324  qaSum+=(*it).second;
2325  if ((*it).second<qaBest) continue;
2326  tkBest=(*it).first; qaBest=(*it).second;
2327  }
2328  dominatrackKey = tkBest; avgQuality = 100*qaBest/(qaSum+1e-10);
2329  for (int hi=0;hi<nHits; hi++) {
2330  const StHit* rHit = recHits[hi];
2331  if (rHit->idTruth()!=tkBest) continue;
2332  DetectorList[rHit->detector()]++;
2333  }
2334  if (DetectorList[kTpcId] > 99) DetectorList[kTpcId] = 99;
2335  if (DetectorList[kSvtId] > 9) DetectorList[kSvtId] = 9;
2336  if (DetectorList[kSsdId] > 9) DetectorList[kSsdId] = 9;
2337  dominatrackHits = DetectorList[kTpcId] + 100*(DetectorList[kSvtId] + 10*DetectorList[kSsdId]);
2338  return;
2339 }
2340 /*
2341  * $Log: StMiniMcMaker.cxx,v $
2342  * Revision 1.50 2019/01/15 19:24:33 genevb
2343  * Kill some memory leaks (thanks, Coverity)
2344  *
2345  * Revision 1.49 2018/01/03 18:18:10 genevb
2346  * idTruths and keys moved from short to int
2347  *
2348  * Revision 1.48 2015/07/29 16:34:31 smirnovd
2349  * Do not store output from function call as it is not used anyway
2350  *
2351  * Revision 1.47 2015/07/29 16:34:24 smirnovd
2352  * Removed defined but not unused local typedefs
2353  *
2354  * Revision 1.46 2015/07/29 16:34:15 smirnovd
2355  * Added std:: to resolve ambiguity for isnan for g++ (4.8)
2356  *
2357  * Revision 1.45 2015/04/30 16:03:31 perev
2358  * Remove redundant automatic CVS comments
2359  *
2360  * Revision 1.44 2015/04/28 16:05:32 perev
2361  * We have changed a default dEdx method in StMuDSTmaker from
2362  * kTruncatedMeanId to kLikelihoodFitId since SL14a.
2363  * With TruncatedMean we are using only 70% of available dE/dx measurements and this is reflected in nHitsDedx
2364  * With LikelihoodFit we are using all available dE/dx measurements.
2365  * Default in StMiniMcEvent is still old (kTruncatedMeanId).
2366  * It is changed now to kLikelihoodFitId and to be back propagated to all releases >= SL14a
2367  * Yuri
2368  *
2369  * Revision 1.43 2014/07/28 17:20:11 jwebb
2370  * Explicit casts from (double) to (float) to satisfy c++ 11 compiler.
2371  *
2372  * Revision 1.42 2013/04/04 21:31:07 perev
2373  * The only MC Pairs added
2374  *
2375  * Revision 1.41 2012/05/25 18:35:30 perev
2376  * Wrong DcaGl fixed. Thanx #2362
2377  *
2378  * Revision 1.40 2012/03/15 23:37:36 perev
2379  * Uncorrected globals added(Chris)
2380  *
2381  * Revision 1.39 2011/10/05 23:07:50 perev
2382  * Commrnt++
2383  *
2384  * Revision 1.38 2011/07/19 19:18:05 perev
2385  * Error handling fixed
2386  *
2387  * Revision 1.37 2011/04/01 20:02:56 perev
2388  * IdTruth part rewritten
2389  *
2390  * Revision 1.36 2011/03/22 00:32:23 perev
2391  * Added impact,phi impact & trigger time
2392  *
2393  * Revision 1.35 2011/02/22 20:42:52 perev
2394  * now int parentParentGeantId
2395  *
2396  * Revision 1.34 2011/02/16 00:50:30 perev
2397  * mPdgId added
2398  *
2399  * Revision 1.33 2010/08/31 20:16:15 fisyak
2400  * Add track seedQuality
2401  *
2402  * Revision 1.32 2010/08/05 15:31:11 jwebb
2403  * Changed the check on valid towers during clustering to suppress the non-
2404  * error issued by StEmcGeom.
2405  *
2406  * Revision 1.31 2010/04/15 19:17:27 fisyak
2407  * Add corrections for AppendMCDaughterTrack from Masayuki Wada
2408  *
2409  * Revision 1.30 2010/01/27 21:28:10 perev
2410  * CleanUp
2411  *
2412  * Revision 1.29 2009/02/02 19:30:50 fisyak
2413  * Set common Hit as no.Tpc + 100*no.Svt + 1000*no.Ssd hits, add protection against empty emcCollection
2414  *
2415  * Revision 1.28 2008/07/17 22:50:52 calderon
2416  * Remove a cut in acceptRaw(const StMcTrack*) that checked on the pseudorapidity.
2417  * This cut was affecting heavy particles thrown flat in rapidity for embedding.
2418  *
2419  * Revision 1.27 2007/12/22 20:31:20 calderon
2420  * Storing of info of 3 EMC towers for each TinyRcTrack and TinyMcTrack.
2421  *
2422  * Revision 1.26 2007/04/26 04:08:05 perev
2423  * Hide StBFChain dependency
2424  *
2425  * Revision 1.25 2007/04/17 05:09:07 perev
2426  * GetTFile()==>StMaker. Jerome request
2427  *
2428  * Revision 1.24 2007/02/23 17:07:41 fisyak
2429  * Resolve bug #682
2430  *
2431  * Revision 1.23 2007/02/21 23:19:36 calderon
2432  * Remove the forcing of the Ghost loop to be "off" by checking the filename.
2433  * This caused problems for embedding when the Ghost loop was requested (the check
2434  * overrode the request). Now, the check only displays a message in the log
2435  * saying that the Ghost loop is turned on.
2436  *
2437  * Revision 1.22 2006/08/05 01:07:45 calderon
2438  * Modified some debugging printouts.
2439  *
2440  * Revision 1.21 2006/07/24 19:04:11 calderon
2441  * Added parent key data member to StTinyMcTrack.
2442  * Added reco key data member to StTinyRcTrack.
2443  * Added code to fill in those data members to StMiniMcMaker. Parent key for
2444  * MC tracks is only filled when track has a valid parent().
2445  *
2446  * Revision 1.20 2006/05/22 18:55:16 calderon
2447  * Changes from the original code by Bum to comply with STAR coding standards.
2448  * First thing is to change the name of the "Helper" file to something that is more in line with the file naming convention.
2449  * This does not fully solve all possible hiccups, because all the functions
2450  * in the "helper" file are defined in global scope.
2451  *
2452  * Revision 1.19 2005/09/29 15:53:08 fisyak
2453  * Persistent StMcEvent
2454  *
2455  * Revision 1.18 2004/07/27 19:34:34 jeromel
2456  * Patch for pileup. Not a primary => a parent but in pileup, not true anymore ...
2457  *
2458  * Revision 1.17 2004/05/03 23:28:39 perev
2459  * double delete of TTree fixed. TFile deletes it himself
2460  *
2461  * Revision 1.16 2004/03/31 23:44:36 calderon
2462  * Function to find the dominatrack, the number of hits belonging to the
2463  * dominatrack and the average hit quality of those hits (based on idTruth and
2464  * quality of StHit).
2465  *
2466  * Revision 1.15 2004/03/30 03:16:15 calderon
2467  * Modifications for running in bfc.
2468  * - Changed to use StiIOInterface (IOMaker in normal mode, TreeMaker in bfc)
2469  * - Cleaned up Init(), InitRun() to handle the changing file names.
2470  * - Initialize lots of variables and pointers in constructor.
2471  * - Delete some pointers in Finish (deleting the TTree causes a seg fault, though.)
2472  * - Note that currently the StHits in the ITTF chain don't have a usedInFit() flag,
2473  * so there will be many messages complaining about this.
2474  * - Removed the mDebug data member, every Maker already has one, so change
2475  * to use that throughout the package.
2476  *
2477  * Revision 1.14 2004/03/15 18:59:47 calderon
2478  * - Added support for encoded common hits. Now the common hits of the TPC and
2479  * the SVT are stored, with the corresponding methods to decode and return these
2480  * values.
2481  * - Added protection for tracks with no particle definition.
2482  * - Added () around call to mcMergedPair[i] to make Insure++ happy.
2483  *
2484  * Revision 1.13 2004/01/26 13:59:26 calderon
2485  * Added the code to fill the global track matches of StMiniMcEvent.
2486  *
2487  * Revision 1.12 2003/09/02 17:58:43 perev
2488  * gcc 3.2 updates + WarnOff
2489  *
2490  * Revision 1.11 2003/07/09 01:07:23 calderon
2491  * Addition of FTPC reference multiplicity
2492  * Addition of other multiplicity values for StMiniMcEvent
2493  * Changes to reflect the use of the setters and getters, no longer
2494  * access the data members directly.
2495  *
2496  * Revision 1.10 2003/05/14 00:12:20 calderon
2497  * The minimc replaces now whatever it finds between the first and last '.', not
2498  * just geant.root, in the creation of the output file name.
2499  *
2500  * Curvature, tan(lambda) and Covariance matrix diagonal elements are now stored
2501  * in rcTrack, for both primary and global tracks.
2502  *
2503  * Revision 1.9 2003/05/08 02:11:43 calderon
2504  * Set the data members for the Svt and Ftpc hits for constMcTrack and for the Svt and Ftpc
2505  * fit points for the RcTrack.
2506  * Use primary track for the fit points in all cases, not the global tracks.
2507  * Selection of West or East Ftpc is based on eta>1.8 or eta<1.8
2508  *
2509  * Revision 1.8 2002/06/28 22:15:12 calderon
2510  * Changes to deal with seg. faults in the file name handling:
2511  * Conventions:
2512  * StMiniMcMaker looks for the input file from the IO maker to figure out
2513  * if the file has changed. This is done using TString::Contains() in Make().
2514  * Usually we will run one file at a time, but in order not to break Bum's scheme of being
2515  * able to process several files in one go, this is left as is. However, for
2516  * embedding, the file name is not enough, in Eric's new scheme there are repeated
2517  * file names. This is resolved by adding a prefix to the output file name. However,
2518  * this prefix should not be overwritten, so the current code only replaces the
2519  * string inside the output file name pertaining to the input file name, and leaves
2520  * the prefix of the output file intact. This was done for embedding looking for
2521  * st_physics, and here is where the problem arose: hijing files begin with a different
2522  * prefix. To solve this problem, the input file name prefix is now an input parameter
2523  * in the macro.
2524  *
2525  * StMiniEmbed.C and StMiniHijing.C now conform to this convention. StMiniEmbed.C
2526  * did not change its prototype, because all embedding files have st_phyics as prefix.
2527  * StMiniHijing.C changed its prototype, now it takes as an input argument the prefix,
2528  * but in order not to break Jenn's scripts if she was already using this macro,
2529  * this parameter was added at the end and defaults to "rcf", which is appropriate
2530  * for hijing files reconstructed in rcf.
2531  *
2532  * Revision 1.7 2002/06/27 17:30:58 jeromel
2533  * Bug fix. NULL+1 caused a crash ...
2534  *
2535  * Revision 1.6 2002/06/11 19:09:35 calderon
2536  * Bug fix: the filename that was set in the macro was being overwritten
2537  * in InitRun, so the emb80x string which was added to the filename was lost.
2538  * This was fixed by not replacing the filename in InitRun and only replacing
2539  * the current filename starting from st_physics.
2540  * and $Id: StMiniMcMaker.cxx,v 1.50 2019/01/15 19:24:33 genevb Exp $ plus header comments for the macros
2541  *
2542  * Revision 1.4 2002/06/06 23:22:34 calderon
2543  * Changes from Jenn:
2544  * -Add needed libs in StMiniHijing.C
2545  * -Properly do an InitRun(int runnumber) method
2546  *
2547 */
Int_t m_Mode
counters
Definition: StMaker.h:81
Bool_t trackOnEmc(StThreeVectorD *position, StThreeVectorD *momentum, const StTrack *const track, Double_t magField, Double_t emcRadius=225.405) const
Track projection utility.
Definition: StHit.h:125
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StMaker.cxx:634
Monte Carlo Track class All information on a simulated track is stored in this class: kinematics...
Definition: StMcTrack.hh:144
double distance(const StThreeVector< double > &p, bool scanPeriods=true) const
minimal distance between point and helix
Definition: StHelix.cc:240
pair< double, double > pathLength(double r) const
path length at given r (cylindrical r)
Definition: StHelix.cc:351
double ycenter() const
x-center of circle in xy-plane
Definition: StHelix.cc:215
Int_t getNextTowerId(Float_t eta, Float_t phi, Int_t nTowersdEta, Int_t nTowersdPhi) const
Return neighbor tower id&#39;s.
virtual Int_t Finish()
Definition: StMaker.cxx:776
rcTpcHitMapType * rcTpcHitMap()
Diff btw global r and phi coords of FTPC hits.
Definition: Stypes.h:44
Filling of StMiniMcEvent classes from StMcEvent, StEvent, StAssociationMaker.
double xcenter() const
aziumth in xy-plane measured from ring center
Definition: StHelix.cc:207
Event data structure to hold all information from a Monte Carlo simulation. This class is the interfa...
Definition: StMcEvent.hh:169
void Clear(Option_t *option="")
User defined functions.
Definition: Stypes.h:41
double Dca(const double point[3], double *dcaErr=0) const
DCA to given space point (with error matrix)