]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetESDFillUnitArrayTracks.cxx
Coding violations corrected.
[u/mrichter/AliRoot.git] / JETAN / AliJetESDFillUnitArrayTracks.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 //======================================================================
18 // ***July 2006
19 // Fill Unit Array class 
20 // Class used by AliJetESDReader to fill a UnitArray from the information extracted 
21 // from the particle tracks
22 // Author: magali.estienne@ires.in2p3.fr
23 //======================================================================
24
25
26 // --- ROOT system ---
27 #include <TVector3.h>
28 #include <TProcessID.h>
29
30 // --- AliRoot header files ---
31 #include "AliJetUnitArray.h"
32 #include "AliJetESDFillUnitArrayTracks.h"
33
34 // --- ROOT system ---
35 class TSystem;
36 class TLorentzVector;
37 class TGeoManager;
38
39 // --- AliRoot header files ---
40 class AliJetFinder;
41 class AliJetReader;
42 class AliJetESDReader;
43
44 ClassImp(AliJetESDFillUnitArrayTracks)
45
46 //_____________________________________________________________________________
47 AliJetESDFillUnitArrayTracks::AliJetESDFillUnitArrayTracks()
48   : AliJetFillUnitArray(),
49     fNumUnits(0),
50     fHadCorr(0),
51     fApplyMIPCorrection(kTRUE),
52     fESD(0x0),
53     fGrid0(0x0),
54     fGrid1(0x0),
55     fGrid2(0x0),
56     fGrid3(0x0),
57     fGrid4(0x0)
58 {
59   // constructor
60 }
61
62 //_____________________________________________________________________________
63 AliJetESDFillUnitArrayTracks::AliJetESDFillUnitArrayTracks(AliESDEvent* esd)
64   : AliJetFillUnitArray(),
65     fNumUnits(0),
66     fHadCorr(0),
67     fApplyMIPCorrection(kTRUE),
68     fESD(esd),
69     fGrid0(0x0),
70     fGrid1(0x0),
71     fGrid2(0x0),
72     fGrid3(0x0),
73     fGrid4(0x0)
74 {
75   // constructor
76 }
77
78 //_____________________________________________________________________________
79 AliJetESDFillUnitArrayTracks::AliJetESDFillUnitArrayTracks(const AliJetESDFillUnitArrayTracks &det)
80   : AliJetFillUnitArray(det),
81     fNumUnits(det.fNumUnits),
82     fHadCorr(det.fHadCorr),
83     fApplyMIPCorrection(det.fApplyMIPCorrection),
84     fESD(det.fESD),
85     fGrid0(det.fGrid0),
86     fGrid1(det.fGrid1),
87     fGrid2(det.fGrid2),
88     fGrid3(det.fGrid3),
89     fGrid4(det.fGrid4)
90 {
91   // Copy constructor
92 }
93
94 //------------------------------------------------------------------
95 AliJetESDFillUnitArrayTracks& AliJetESDFillUnitArrayTracks::operator=(const AliJetESDFillUnitArrayTracks& other)
96 {
97   // Assignment
98
99     fNumUnits = other.fNumUnits;
100     fHadCorr = other.fHadCorr;
101     fApplyMIPCorrection = other.fApplyMIPCorrection;
102     fESD = other.fESD;
103     fGrid0 = other.fGrid0;
104     fGrid1 = other.fGrid1;
105     fGrid2 = other.fGrid2;
106     fGrid3 = other.fGrid3;
107     fGrid4 = other.fGrid4;
108
109     return (*this);
110 }
111
112 //____________________________________________________________________________
113 void AliJetESDFillUnitArrayTracks::InitParameters()
114 {
115   fHadCorr        = 0;                 // For hadron correction
116   fNumUnits = fGeom->GetNCells();      // Number of towers in EMCAL
117   cout << "In AliJetESDFillUnitArrayTracks:InitParameters(), Ncells : " << fNumUnits << endl;
118
119   fTPCGrid->GetAccParam(fNphi,fNeta,fPhiMin, 
120                         fPhiMax,fEtaMin,fEtaMax);
121   fTPCGrid->GetBinParam(fPhiBinInTPCAcc,fEtaBinInTPCAcc, 
122                         fPhiBinInEMCalAcc,fEtaBinInEMCalAcc,fNbinPhi);
123
124   fIndex = fTPCGrid->GetIndexObject();
125
126   if(fDebug>20){
127     for(Int_t i=0; i<fNphi+1; i++)
128       for(Int_t j=0; j<fNeta+1; j++) {cout << "fIndex[" << i << "," << j << "] : " <<
129           (*fIndex)(i,j) << endl; }
130   } 
131   if(fDebug>1) printf("\n Parameters initiated ! \n");
132 }
133
134 //_____________________________________________________________________________
135 AliJetESDFillUnitArrayTracks::~AliJetESDFillUnitArrayTracks()
136 {
137   // destructor
138 }
139
140 //_____________________________________________________________________________
141 void AliJetESDFillUnitArrayTracks::Exec(Option_t* const /*option*/)
142 {
143   //
144   // Main method.
145   // Fill the unit array with the charged particle information in ESD
146   //
147
148   fDebug = fReaderHeader->GetDebug();
149   
150   // Set parameters
151   InitParameters();
152
153   // get number of tracks in event (for the loop)
154   Int_t goodTrack = 0;
155   Int_t nt = 0;
156 //(not used ?)  Int_t nt2 = 0;
157   Int_t nmax = 0;
158   Float_t pt,eta,phi;
159   //  Int_t sflag = 0;
160   TVector3 p3;
161
162   if(fDebug>1) cout << "fESD in Fill array : " << fESD << endl;
163
164   nt = fESD->GetNumberOfTracks();
165   if(fDebug>1) cout << "Number of Tracks in ESD : " << nt << endl;
166  
167   // temporary storage of signal and pt cut flag
168   Int_t* sflag  = new Int_t[nt];
169   Int_t* cflag  = new Int_t[nt];
170
171   // get cuts set by user
172   Float_t ptMin  = fReaderHeader->GetPtCut();
173   Float_t etaMin = fReaderHeader->GetFiducialEtaMin();
174   Float_t etaMax = fReaderHeader->GetFiducialEtaMax();  
175   fOpt = fReaderHeader->GetDetector();
176   fDZ  = fReaderHeader->GetDZ();
177
178   Int_t nTracksEmcal      = 0;
179   Int_t nTracksEmcalDZ    = 0;
180   Int_t nTracksTpc        = 0;
181   Int_t nTracksTpcOnly    = 0;
182   Int_t nTracksEmcalCut   = 0;
183   Int_t nTracksEmcalDZCut = 0;
184   Int_t nTracksTpcCut     = 0;
185   Int_t nTracksTpcOnlyCut = 0;
186
187 //(not used ?)  Int_t nElements  = fTPCGrid->GetNEntries();
188 //(not used ?)  Int_t nElements2 = fEMCalGrid->GetNEntries();
189   fGrid = fTPCGrid->GetGridType();
190
191 //(not used ?)  Int_t nEntries = fUnitArray->GetEntries();
192
193 //(not used ?)  Int_t nRefEnt = fRefArray->GetEntries();
194
195   //loop over tracks
196     nmax = nt;  
197     for (Int_t it = 0; it < nmax; it++) {
198       AliESDtrack *track;
199       track = fESD->GetTrack(it);
200       UInt_t status = track->GetStatus();
201     
202       Double_t mom[3];
203       track->GetPxPyPz(mom);
204
205       p3.SetXYZ(mom[0],mom[1],mom[2]);
206       pt = p3.Pt();
207
208       Float_t mass = 0.;
209       mass = track->GetMass();
210
211       if ((status & AliESDtrack::kTPCrefit) == 0)    continue;      // quality check
212       if ((status & AliESDtrack::kITSrefit) == 0)    continue;      // quality check
213       if (((AliJetESDReaderHeader*) fReaderHeader)->ReadSignalOnly() 
214           && TMath::Abs(track->GetLabel()) > 10000)  continue;   // quality check
215       if (((AliJetESDReaderHeader*) fReaderHeader)->ReadBkgdOnly() 
216           && TMath::Abs(track->GetLabel()) < 10000)  continue;   // quality check
217        eta = p3.Eta();
218        phi = ( (p3.Phi()) < 0) ? (p3.Phi()) + 2. * TMath::Pi() : (p3.Phi());
219    
220       if ( (eta > etaMax) || (eta < etaMin)) continue;           // checking eta cut
221
222       sflag[goodTrack]=0;
223       if (TMath::Abs(track->GetLabel()) < 10000) sflag[goodTrack]=1;
224       cflag[goodTrack]=0;
225       if (pt > ptMin) cflag[goodTrack]=1;                       // pt cut
226
227       if(fGrid==0)
228         {
229           // Only TPC filled from its grid in its total acceptance
230           
231           Int_t idTPC = fTPCGrid->GetIndex(phi,eta);
232           Bool_t ok = kFALSE;
233           Bool_t ref = kFALSE;
234
235           AliJetUnitArray *uArray = (AliJetUnitArray*)fUnitArray->At(idTPC-1);
236           TRefArray *reference = uArray->GetUnitTrackRef();
237           reference->Add(track);
238
239           Float_t unitEnergy = 0.;
240           unitEnergy = uArray->GetUnitEnergy();
241           // nTracksTpcOnly is necessary to count the number of candidate cells
242           // but it doesn't give the real multiplicity -> it will be extracted 
243           // in the jet finder through the reference to tracks
244           if(unitEnergy==0.){
245             nTracksTpcOnly++;
246             ok = kTRUE;
247             ref = kTRUE;
248           }
249
250           // Fill energy in TPC acceptance
251           uArray->SetUnitEnergy(unitEnergy + pt);
252           uArray->SetUnitMass(mass);
253
254           // Pt cut flag
255           if(uArray->GetUnitEnergy()<ptMin){
256             uArray->SetUnitCutFlag(kPtSmaller);
257           }
258           else {
259             uArray->SetUnitCutFlag(kPtHigher);
260             if(ok) nTracksTpcOnlyCut++;
261           }
262
263           // Signal flag
264           if(sflag[goodTrack] == 1) { 
265             uArray->SetUnitSignalFlag(kGood);
266             uArray->SetUnitSignalFlagC(kFALSE,kGood);
267           } else uArray->SetUnitSignalFlagC(kFALSE,kBad);
268
269           if(uArray->GetUnitEnergy()>0 && ref){
270             if(!fProcId) {
271               fProcId = kTRUE;
272               new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray));
273             }
274             fRefArray->Add(uArray);
275           }
276         }
277     
278       if(fGrid==1)
279         {
280           Int_t nElements = fTPCGrid->GetNEntries();
281           // Fill track information in EMCAL acceptance
282             if((eta >= fEtaMin && eta <= fEtaMax) &&
283              (phi >= fPhiMin && phi <= fPhiMax))// &&
284             {
285
286               // Include dead-zones
287               if(fDZ)
288                 {
289                   Double_t phimin0 = 0., phimin1 = 0., phimin2 = 0., phimin3 = 0., phimin4 = 0.;
290                   Double_t phimax0 = 0., phimax1 = 0., phimax2 = 0., phimax3 = 0., phimax4 = 0.;
291                   fGeom->GetPhiBoundariesOfSMGap(0,phimin0,phimax0);
292                   fGeom->GetPhiBoundariesOfSMGap(1,phimin1,phimax1);
293                   fGeom->GetPhiBoundariesOfSMGap(2,phimin2,phimax2);
294                   fGeom->GetPhiBoundariesOfSMGap(3,phimin3,phimax3);
295                   fGeom->GetPhiBoundariesOfSMGap(4,phimin4,phimax4);
296                   Int_t n0 = fGrid0->GetNEntries();
297                   Int_t n1 = fGrid1->GetNEntries();
298                   Int_t n2 = fGrid2->GetNEntries();
299                   Int_t n3 = fGrid3->GetNEntries();
300 //(not used ?)            Int_t n4 = fGrid4->GetNEntries();
301                   
302                   if(phi >= phimin0 && phi <= phimax0){
303                     Int_t id0 = fGrid0->GetIndex(phi,eta)-1;
304                     AliJetUnitArray *uArray0 = (AliJetUnitArray*)fUnitArray->At(id0+fNumUnits+nElements);
305                     TRefArray *reference0 = uArray0->GetUnitTrackRef();
306                     reference0->Add(track);
307
308                     Float_t uEnergy0 = uArray0->GetUnitEnergy();
309                     Bool_t ok0 = kFALSE;
310                     Bool_t ref0 = kFALSE;
311                     if(uEnergy0==0.){
312                       nTracksEmcalDZ++;
313                       ok0 = kTRUE;
314                       ref0 = kTRUE;
315                     }
316                     uArray0->SetUnitEnergy(uEnergy0+pt);
317
318                     if(uArray0->GetUnitEnergy()<ptMin)
319                       uArray0->SetUnitCutFlag(kPtSmaller);
320                     else {
321                       uArray0->SetUnitCutFlag(kPtHigher);
322                       if(ok0) nTracksEmcalDZCut++;
323                     }
324                     if(sflag[goodTrack] == 1) {
325                       uArray0->SetUnitSignalFlag(kGood);
326                       uArray0->SetUnitSignalFlagC(kFALSE,kGood);
327                     } else uArray0->SetUnitSignalFlagC(kFALSE,kBad);
328
329                     if(uArray0->GetUnitEnergy()>0 && ref0){
330                       if(!fProcId){
331                         new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray0));
332                         fProcId = kTRUE;
333                       }
334                       fRefArray->Add(uArray0);
335                     }
336                   }
337
338                   if(phi >= phimin1 && phi <= phimax1){
339                     Int_t id1 = fGrid1->GetIndex(phi,eta)-1+n0;
340                     AliJetUnitArray *uArray1 = (AliJetUnitArray*)fUnitArray->At(id1+fNumUnits+nElements);
341                     TRefArray *reference1 = uArray1->GetUnitTrackRef();
342                     reference1->Add(track);
343
344                     Float_t uEnergy1 = uArray1->GetUnitEnergy();
345                     Bool_t ok1 = kFALSE;
346                     Bool_t ref1 = kFALSE;
347                     if(uEnergy1==0.){
348                       nTracksEmcalDZ++;
349                       ok1 = kTRUE;
350                       ref1 = kTRUE;
351                     }
352                     uArray1->SetUnitEnergy(uEnergy1+pt);
353
354                     if(uArray1->GetUnitEnergy()<ptMin)
355                       uArray1->SetUnitCutFlag(kPtSmaller);
356                     else {
357                       uArray1->SetUnitCutFlag(kPtHigher);
358                       if(ok1) nTracksEmcalDZCut++;
359                     }
360                     if(sflag[goodTrack] == 1) {
361                       uArray1->SetUnitSignalFlag(kGood);
362                       uArray1->SetUnitSignalFlagC(kFALSE,kGood);
363                     } else uArray1->SetUnitSignalFlagC(kFALSE,kBad);
364
365                     if(uArray1->GetUnitEnergy()>0 && ref1){
366                       if(!fProcId){
367                         new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray1));
368                         fProcId = kTRUE;
369                       }
370                       fRefArray->Add(uArray1);
371                     }
372                   }
373
374                   if(phi >= phimin2 && phi <= phimax2){
375                     Int_t id2 = fGrid2->GetIndex(phi,eta)-1+n0+n1;
376                     AliJetUnitArray *uArray2 = (AliJetUnitArray*)fUnitArray->At(id2+fNumUnits+nElements);
377                     TRefArray *reference2 = uArray2->GetUnitTrackRef();
378                     reference2->Add(track);
379
380                     Float_t uEnergy2 = uArray2->GetUnitEnergy();
381                     Bool_t ok2 = kFALSE;
382                     Bool_t ref2 = kFALSE;
383                     if(uEnergy2==0.){
384                       nTracksEmcalDZ++;
385                       ok2 = kTRUE;
386                       ref2 = kTRUE;
387                     }
388                     uArray2->SetUnitEnergy(uEnergy2+pt);
389
390                     if(uArray2->GetUnitEnergy()<ptMin)
391                       uArray2->SetUnitCutFlag(kPtSmaller);
392                     else {
393                       uArray2->SetUnitCutFlag(kPtHigher);
394                       if(ok2) nTracksEmcalDZCut++;
395                     }
396                     if(sflag[goodTrack] == 1) {
397                       uArray2->SetUnitSignalFlag(kGood);
398                       uArray2->SetUnitSignalFlagC(kFALSE,kGood);
399                     } else uArray2->SetUnitSignalFlagC(kFALSE,kBad);
400
401                     if(uArray2->GetUnitEnergy()>0 && ref2){
402                       if(!fProcId){
403                         new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray2));
404                         fProcId = kTRUE;
405                       }
406                       fRefArray->Add(uArray2);
407                     }
408                   }
409
410                   if(phi >= phimin3 && phi <= phimax3){
411                     Int_t id3 = fGrid3->GetIndex(phi,eta)-1+n0+n1+n2;
412                     AliJetUnitArray *uArray3 = (AliJetUnitArray*)fUnitArray->At(id3+fNumUnits+nElements);
413                     TRefArray *reference3 = uArray3->GetUnitTrackRef();
414                     reference3->Add(track);
415
416                     Float_t uEnergy3 = uArray3->GetUnitEnergy();
417                     Bool_t ok3 = kFALSE;
418                     Bool_t ref3 = kFALSE;
419                     if(uEnergy3==0.){
420                       nTracksEmcalDZ++;
421                       ok3 = kTRUE;
422                       ref3 = kTRUE;
423                     }
424                     uArray3->SetUnitEnergy(uEnergy3+pt);
425
426                     if(uArray3->GetUnitEnergy()<ptMin)
427                       uArray3->SetUnitCutFlag(kPtSmaller);
428                     else {
429                       uArray3->SetUnitCutFlag(kPtHigher);
430                       if(ok3) nTracksEmcalDZCut++;
431                     }
432                     if(sflag[goodTrack] == 1) {
433                       uArray3->SetUnitSignalFlag(kGood);
434                       uArray3->SetUnitSignalFlagC(kFALSE,kGood);
435                     } else uArray3->SetUnitSignalFlagC(kFALSE,kBad);
436
437                     if(uArray3->GetUnitEnergy()>0 && ref3){
438                       if(!fProcId){
439                         new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray3));
440                         fProcId = kTRUE;
441                       }
442                       fRefArray->Add(uArray3);
443                     }
444                   }
445
446                   if(phi >= phimin4 && phi <= phimax4){
447                     Int_t id4 = fGrid4->GetIndex(phi,eta)-1+n0+n1+n2+n3;
448                     AliJetUnitArray *uArray4 = (AliJetUnitArray*)fUnitArray->At(id4+fNumUnits+nElements);
449                     TRefArray *reference4 = uArray4->GetUnitTrackRef();
450                     reference4->Add(track);
451
452                     Float_t uEnergy4 = uArray4->GetUnitEnergy();
453                     Bool_t ok4 = kFALSE;
454                     Bool_t ref4 = kFALSE;
455                     if(uEnergy4==0.){
456                       nTracksEmcalDZ++;
457                       ok4 = kTRUE;
458                       ref4 = kTRUE;
459                     }
460                     uArray4->SetUnitEnergy(uEnergy4+pt);
461
462                     if(uArray4->GetUnitEnergy()<ptMin)
463                       uArray4->SetUnitCutFlag(kPtSmaller);
464                     else {
465                       uArray4->SetUnitCutFlag(kPtHigher);
466                       if(ok4) nTracksEmcalDZCut++;
467                     }
468                     if(sflag[goodTrack] == 1) {
469                       uArray4->SetUnitSignalFlag(kGood);
470                       uArray4->SetUnitSignalFlagC(kFALSE,kGood);
471                     } else uArray4->SetUnitSignalFlagC(kFALSE,kBad);
472
473                     if(uArray4->GetUnitEnergy()>0 && ref4){
474                       if(!fProcId){
475                         new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray4));
476                         fProcId = kTRUE;
477                       }
478                       fRefArray->Add(uArray4);
479                     }
480                   }
481                 } // end fDZ
482             
483               Int_t towerID = 0;
484               fGeom->GetAbsCellIdFromEtaPhi(eta,phi,towerID);
485               if(towerID==-1) continue;
486               
487               AliJetUnitArray *uArray = (AliJetUnitArray*)fUnitArray->At(towerID);
488               TRefArray *reference = uArray->GetUnitTrackRef();
489               reference->Add(track);
490
491               Float_t unitEnergy = uArray->GetUnitEnergy(); 
492               Bool_t ok = kFALSE;
493               Bool_t ref = kFALSE;
494               if(unitEnergy==0.){
495                 nTracksEmcal++;
496                 ok=kTRUE;
497                 ref=kTRUE;
498               }
499
500               // Do Hadron Correction
501               // This is under construction !!!!!!!!!!!!!!!!!!!!!!!
502               // Parametrization to be added
503               if (fApplyMIPCorrection != 0) 
504                 { 
505 //                Float_t   hCEnergy = fHadCorr->GetEnergy(p3.Mag(), (Double_t)eta,0);
506 //                unitEnergy -= hCEnergy*TMath::Sin(2.0*TMath::ATan(TMath::Exp(-eta)));
507                 } //end Hadron Correction loop
508
509               uArray->SetUnitEnergy(unitEnergy + pt);
510
511               // Put a pt cut flag
512               if(uArray->GetUnitEnergy()<ptMin){
513                 uArray->SetUnitCutFlag(kPtSmaller);
514               }
515               else {
516                 uArray->SetUnitCutFlag(kPtHigher);
517                 if(ok) nTracksEmcalCut++;
518               }
519
520               // Signal flag
521               if(sflag[goodTrack] == 1) { 
522                 uArray->SetUnitSignalFlag(kGood);
523                 uArray->SetUnitSignalFlagC(kFALSE,kGood);
524               } else uArray->SetUnitSignalFlagC(kFALSE,kBad);
525
526             
527               if(uArray->GetUnitEnergy()>0 && ref){
528                 if(!fProcId){
529                   new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray));
530                   fProcId = kTRUE;
531                 }
532                 fRefArray->Add(uArray);
533               }
534             } // end loop on EMCal acceptance cut + tracks quality
535           else{ 
536             // Outside EMCal acceptance
537           
538             Int_t idTPC = fTPCGrid->GetIndex(phi,eta);
539
540             AliJetUnitArray *uArray = (AliJetUnitArray*)fUnitArray->At(fNumUnits-1+idTPC);
541             TRefArray *reference = uArray->GetUnitTrackRef();
542             reference->Add(track);
543
544             Float_t unitEnergy2 = uArray->GetUnitEnergy(); 
545             Bool_t okout = kFALSE;
546             Bool_t refout = kFALSE;
547             if(unitEnergy2==0.){
548               nTracksTpc++;
549               okout=kTRUE;
550               refout=kTRUE;
551             }
552             // Fill energy outside emcal acceptance
553             uArray->SetUnitEnergy(unitEnergy2 + pt);
554           
555             // Pt cut flag
556             if(uArray->GetUnitEnergy()<ptMin){
557               uArray->SetUnitCutFlag(kPtSmaller);
558             }
559             else {
560               uArray->SetUnitCutFlag(kPtHigher);
561               if(okout) nTracksTpcCut++;
562             }
563             // Signal flag
564             if(sflag[goodTrack] == 1) {
565               uArray->SetUnitSignalFlag(kGood);
566               uArray->SetUnitSignalFlagC(kFALSE,kGood);
567             } else uArray->SetUnitSignalFlagC(kFALSE,kBad);
568
569
570             if(uArray->GetUnitEnergy()>0 && refout){
571               if(!fProcId){
572                 new(fRefArray) TRefArray(TProcessID::GetProcessWithUID(uArray));
573                 fProcId = kTRUE;
574               }
575               fRefArray->Add(uArray);
576             }
577           }
578         } // end fGrid==1
579
580       goodTrack++;
581
582     } // end loop on entries (tpc tracks)
583
584   if(fDebug>0)
585     {
586       cout << "End of Tracks %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%" << endl;
587       cout << "goodTracks: " << goodTrack << endl;
588     }
589
590   delete[] sflag;
591   delete[] cflag;
592       
593   if(fGrid==0) {
594     fNTracks = nTracksTpcOnly;
595     fNTracksCut = nTracksTpcOnlyCut;
596     if(fDebug>10){
597       cout << "fNTracks : " << fNTracks << endl;
598       cout << "fNTracksCut : " << fNTracksCut << endl;
599     }
600   }
601   if(fGrid==1) {
602     fNTracks = nTracksEmcal+nTracksEmcalDZ+nTracksTpc;
603     fNTracksCut = nTracksEmcalCut+nTracksEmcalDZCut+nTracksTpcCut;
604     if(fDebug>10){
605       cout << "fNTracks : " << fNTracks << endl;
606       cout << "fNTracksCut : " << fNTracksCut << endl;
607     }
608   }  
609   
610 }
611
612
613
614
615
616
617
618
619
620
621
622
623
624