]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHTracker.cxx
ca731514b316c44c80beeae26e7cd18bf4d5c568
[u/mrichter/AliRoot.git] / RICH / AliRICHTracker.cxx
1 #include "AliRICHTracker.h" //class header
2 #include "AliRICH.h"
3 #include "AliRICHRecon.h"
4 #include <AliESD.h>
5 #include <TVector3.h>
6 #include <TTree.h>          //EsdPrint() 
7 #include <TFile.h>          //EsdPrint()    
8 #include "AliRICHHelix.h"
9 #include <AliMagF.h>
10 #include <AliStack.h>
11 #include <TParticle.h>
12 #include <TMath.h>
13 #include <AliRun.h>
14 #include <TNtupleD.h>            //RecWithStack();
15 #include <AliTrackPointArray.h>  //GetTrackPoint()
16 #include <AliAlignObj.h>         //GetTrackPoint()
17 ClassImp(AliRICHTracker)
18 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
19 AliRICHTracker::AliRICHTracker():AliTracker()
20 {
21 // AliRICHTracker is created from AliReconstraction::Run() which invokes AliReconstraction::CreateTrackers() 
22 // which in turn invokes AliRICHReconstructor::CreateTracker(). 
23 // Note that this is done just once per session before AliReconstruction::Run() goes to events loop.
24   AliRICHParam::Instance()->CdbRead(0,0); 
25   for(Int_t i=0;i<5;i++)fErrPar[i]=0;
26 }
27 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
28 Bool_t AliRICHTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
29 {
30 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
31 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
32 // Arguments: idx- cluster index which is stored by RICH in AliESDtrack
33 //            point- reference to the object where to store the point     
34 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
35   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
36   Int_t iCham=idx/1000000;
37   Int_t iClu=idx%1000000;
38   point.SetVolumeID(AliAlignObj::LayerToVolUID(AliAlignObj::kRICH,iCham-1));//layer and chamber number
39   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));  
40   AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iCham)->UncheckedAt(iClu);//get pointer to cluster
41   TVector3 mars=AliRICHParam::Instance()->Lors2Mars(iCham,pClu->X(),pClu->Y());
42   point.SetXYZ(mars.X(),mars.Y(),mars.Z());
43   return kTRUE;
44 }
45 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
46 Int_t AliRICHTracker::LoadClusters(TTree *pCluTree)
47 {
48 // Interface callback methode invoked from AliReconstruction::RunTracking() to load RICH clusters for RICH
49 // Arguments: pCluTree- pointer to clusters tree got by AliRICHLoader::LoadRecPoints("read") then AliRICHLoader::TreeR()
50 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
51   AliDebug(1,"Start.");  pCluTree->GetEntry(0);  AliDebug(1,"Stop."); return 0;
52 }
53 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
54 Int_t AliRICHTracker::PropagateBack(AliESD *pESD)
55 {
56 // Interface callback methode invoked by AliRecontruction::RunTracking() during tracking after TOF. It's done just once per event
57 // Arguments: pESD - pointer to Event Summary Data class instance which contains a list of tracks
58 //   Returns: error code, 0 if no errors   
59   Int_t iNtracks=pESD->GetNumberOfTracks();
60   AliDebug(1,Form("Start with %i tracks",iNtracks));
61   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));  
62   AliRICHRecon recon;
63   AliPID pid; // needed to retrive all the PID info
64    
65   for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
66     AliESDtrack *pTrack = pESD->GetTrack(iTrk);// get next reconstructed track    
67     AliRICHHelix helix(pTrack->X3(),pTrack->P3(),(Int_t)pTrack->GetSign(),-0.1*GetBz()); //construct helix out of track running parameters  
68      //Printf(" magnetic field %f charged %f\n",GetBz(),pTrack->GetSign()); helix.Print("Track");
69     Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());    
70     if(!iChamber) continue;                                         //no intersection with chambers, ignore this track go after the next one
71   
72     //find MIP cluster candidate (closest to track intersection point cluster with large enough QDC)
73     Double_t    dR=9999,   dX=9999,   dY=9999;                     //distance between track-PC intersection point and current cluster
74     Double_t mipDr=9999,mipDx=9999,mipDy=9999,mipX=9999,mipY=9999; //nearest cluster parameters
75     Int_t   iMipId=-1;                         //index of this nearest cluster
76     for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
77       AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
78       if(pClu->Q()<AliRICHParam::QthMIP()) continue; //to low QDC, go after another one
79       pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY); //get distance for current cluster
80       if(dR<mipDr){iMipId=iClu; mipDr=dR; mipDx=dX; mipDy=dY; mipX=pClu->X(); mipY=pClu->Y();} //current cluster is closer, overwrite data for min cluster
81     }//clusters loop for intersected chamber
82
83     pTrack->SetRICHthetaPhi(helix.Ploc().Theta(),helix.Ploc().Phi()); //store track impact angles with respect to RICH planes
84     pTrack->SetRICHdxdy(mipDx,mipDy);                                 //distance between track-PC intersection and closest cluster with Qdc>100
85     pTrack->SetRICHmipXY(mipX,mipY);                                  //position of that closest cluster with Qdc>100
86     
87     if(iMipId==-1)                        {pTrack->SetRICHsignal(kMipQdcCut);  continue;} //no cluster with enough QDC found
88     if(mipDr>AliRICHParam::DmatchMIP())   {pTrack->SetRICHsignal(kMipDistCut); continue;} //closest cluster with enough carge is still too far 
89   
90     pTrack->SetRICHcluster(iMipId+1000000*iChamber);                                //set mip cluster index
91     pTrack->SetRICHsignal(recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId));//search for mean Cerenkov angle for this track
92     pTrack->SetRICHnclusters(iMipId);                                               //on return iMipId is number of photon clusters accepted in reconstruction
93
94     AliDebug(1,Form("Ch=%i PC Intersection=(%5.2f,%5.2f) cm MIP cluster dist=(%5.2f,%5.2f)=%5.2f cm ThetaCkov=%f",
95                      iChamber,helix.PosPc().X(),helix.PosPc().Y(),            mipDx,mipDy,mipDr,     pTrack->GetRICHsignal()));
96     
97 //here comes PID calculations    
98     if(pTrack->GetRICHsignal()>0) {
99       AliDebug(1,Form("Start to assign the probabilities"));
100       Double_t sigmaPID[AliPID::kSPECIES];  Double_t richPID[AliPID::kSPECIES];
101       for (Int_t iPid=0;iPid<AliPID::kSPECIES;iPid++){//PID loop
102         sigmaPID[iPid] = 0; fErrPar[iPid] = 0;
103         for(Int_t iCkov=0;iCkov<pRich->Clus(iChamber)->GetEntries();iCkov++){//Ckov candidates loop ????????????? somehting fomr AliRICHRecon must be
104           recon.SetPhotonIndex(iCkov);
105           if(recon.GetPhotonFlag() == 2) {
106             Double_t thetaCkov=0.6; //??????????????????
107             Double_t   phiCkov=0.6; //??????????????????
108             Double_t  thetaTrk=recon.GetTrackTheta();
109             Double_t    phiTrk=(recon.GetPhiPoint()-recon.GetTrackPhi());
110             Double_t     beta =pTrack->GetP()/TMath::Sqrt((pTrack->GetP()*pTrack->GetP()+pid.ParticleMass(iPid)*pid.ParticleMass(iPid)));
111             Double_t   sigma2 = AliRICHParam::SigmaSinglePhotonFormula(thetaCkov,phiCkov,thetaTrk,phiTrk,beta);
112             if(sigma2>0) sigmaPID[iPid] += 1/sigma2;
113           }
114         }//Ckov candidates loop
115         
116         if (sigmaPID[iPid]>0)
117           sigmaPID[iPid] *= (Double_t)(iMipId-recon.GetPhotBKG())/(Double_t)(iMipId); // n total phots, m are background...the sigma are scaled..
118           if(sigmaPID[iPid]>0) sigmaPID[iPid] = 1/TMath::Sqrt(sigmaPID[iPid])*0.001; // sigma from parametrization are in mrad...
119           else                  sigmaPID[iPid] = 0;
120           fErrPar[iPid]=sigmaPID[iPid];
121         AliDebug(1,Form("sigma for %s is %f rad",AliPID::ParticleName(iPid),sigmaPID[iPid]));
122       }
123       CalcProb(pTrack->GetRICHsignal(),pTrack->GetP(),sigmaPID,richPID);
124       pTrack->SetRICHpid(richPID);         
125       AliDebug(1,Form("PROBABILITIES ---> %f - %f - %f - %f - %f",richPID[0],richPID[1],richPID[2],richPID[3],richPID[4]));
126     }//if(pTrack->GetRICHsignal())    
127   }//ESD tracks loop
128   AliDebug(1,"Stop pattern recognition");
129   return 0; // error code: 0=no error;
130 }//PropagateBack()
131 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
132 void AliRICHTracker::RecWithStack(TNtupleD *hn)
133 {
134 // Reconstruction for particles from STACK. This methode is to be used for RICH standalone when no other detectors are switched on, so normal tracking is not available.
135 // Arguments: hn- output ntuple where to store all variables
136 //   Returns: none       
137   AliDebug(1,"Start.");  
138   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
139   
140 //  pRich->GetLoader()->GetRunLoader()->LoadHeader();
141   if(!pRich->GetLoader()->GetRunLoader()->TreeK()) pRich->GetLoader()->GetRunLoader()->LoadKinematics();
142   AliStack *pStack =   pRich->GetLoader()->GetRunLoader()->Stack();
143   if(!pStack) {AliDebug(1,Form("No STACK found in AliRoot"));return;}
144   Int_t iNtracks=pStack->GetNtrack();
145   AliDebug(1,Form(" Start reconstruction with %i track(s) from Stack",iNtracks));
146   
147   Double_t hnvec[20];
148   
149   Double_t b=GetFieldMap()->SolenoidField()/10;// magnetic field in Tesla
150   AliDebug(1,Form("Start with simulated %i tracks in %f Tesla field",iNtracks,b));
151   TVector3 x0(0,0,0); TVector3 p0(0,0,0);//tmp storage for AliRICHHelix
152   
153
154   if(pRich->GetLoader()->LoadRecPoints()) {AliDebug(1,Form("No clusters found in RICH"));return;}
155   pRich->GetLoader()->TreeR()->GetEntry(0);
156
157   AliRICHRecon recon;
158   for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//stack particles loop
159     TParticle *pParticle = pStack->Particle(iTrackN);
160     if(!pParticle) {AliDebug(1,Form("Not a valid TParticle pointer. Track skipped"));continue;}
161     AliDebug(1,Form(" PDG code : %i",pParticle->GetPdgCode()));
162 //
163 // problem of PDG code of some extra particles to be solved!!!!!!!!!
164 //
165 // found problem! Look in TRD directory : codes from Fluka are :
166 //
167 //    if ((pdg_code == 10010020) ||
168 //        (pdg_code == 10010030) ||
169 //        (pdg_code == 50000050) ||
170 //        (pdg_code == 50000051) ||
171 //        (pdg_code == 10020040)) {
172 //
173     if(pParticle->GetPdgCode()>=50000050||pParticle->GetPdgCode()==0||pParticle->GetPdgCode()>10000) {AliDebug(1,Form("A photon as track... Track skipped"));continue;}
174     if(!pParticle->GetPDG()) continue;
175 //
176 // to be updated for us!!
177 //
178     AliDebug(1,Form("Track %i is a %s with charge %i and momentum %f",
179             iTrackN,pParticle->GetPDG()->GetName(),Int_t(pParticle->GetPDG()->Charge()),pParticle->P()));
180 //    if(pParticle->GetMother(0)!=-1) continue; //consider only primaries
181     if(pParticle->GetPDG()->Charge()==0||TMath::Abs(Int_t(pParticle->GetPDG()->Charge()))!=3) continue; //to avoid photons from stack...
182     hnvec[0]=pParticle->P();
183     hnvec[1]=pParticle->GetPDG()->Charge();
184     
185     p0.SetMagThetaPhi(pParticle->P(),pParticle->Theta(),pParticle->Phi());   x0.SetXYZ(pParticle->Vx(),pParticle->Vy(),pParticle->Vz());
186     AliRICHHelix helix(x0,p0,TMath::Sign(1,(Int_t)pParticle->GetPDG()->Charge()),b);   
187     Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());        
188     if(!iChamber) continue;// no intersection with RICH found
189     
190     hnvec[2]=helix.Ploc().Theta();
191     hnvec[3]=helix.Ploc().Phi();
192     hnvec[4]=helix.PosPc().X();
193     hnvec[5]=helix.PosPc().Y();
194     
195     Double_t dX,dY,dR,dRmip=9999;      //min distance between clusters and track position on PC 
196     Int_t iMipId=-1;       //index of that min distance cluster 
197     for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
198       AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
199       pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY);//ditance between current cluster and helix intersection with PC
200       if(dR<dRmip){dRmip=dR; hnvec[6]=pClu->X();hnvec[7]=pClu->Y();hnvec[8]=pClu->Q();
201                                                   iMipId=1000000*iChamber+iClu;}//find cluster nearest to the track 
202       
203     }//clusters loop for intersected chamber
204     
205     
206     hnvec[9]=recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId); //search for mean Cerenkov angle for this track
207     hnvec[10]=iMipId;//on return from ThetaCerenkov() contains number of photon candidates accepted
208     hnvec[11]=(Double_t)iMipId;
209     hnvec[12]=(Double_t)iChamber;
210     hnvec[13]=(Double_t)pParticle->GetPdgCode();
211     if(hn) hn->Fill(hnvec);
212     AliDebug(1,Form("FINAL Theta Cerenkov=%f",hnvec[9]));
213   }//stack particles loop
214   
215   pRich->GetLoader()->UnloadRecPoints();
216   AliDebug(1,"Stop.");  
217 }//RecWithStack
218 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
219 void AliRICHTracker::CalcProb(Double_t thetaCer,Double_t pmod, Double_t *sigmaPID, Double_t *richPID)
220 {
221 // Calculates probability to be a electron-muon-pion-kaon-proton
222 // from the given Cerenkov angle and momentum assuming no initial particle composition
223 // (i.e. apriory probability to be the particle of the given sort is the same for all sorts)
224   Double_t height[AliPID::kSPECIES];Double_t totalHeight=0;
225   Double_t thetaTh[AliPID::kSPECIES];
226   for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
227     height[iPart]=0;
228     Double_t mass = AliRICHParam::fgMass[iPart];
229     Double_t refIndex=AliRICHParam::Instance()->IdxC6F14(AliRICHParam::EckovMean());
230     Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(refIndex*pmod);
231     thetaTh[iPart]=0;
232     if(cosThetaTh>=1) continue;
233     thetaTh[iPart] = TMath::ACos(cosThetaTh);
234 //    Double_t sinThetaThNorm = TMath::Sin(thetaTh)/TMath::Sqrt(1-1/(refIndex*refIndex));
235 //    Double_t sigmaThetaTh = (0.014*(1/sinThetaThNorm-1) + 0.0043)*1.25;
236 //    height[iPart] = TMath::Gaus(thetaCer,thetaTh,sigmaThetaTh);
237     if(sigmaPID[iPart]>0) height[iPart] = TMath::Gaus(thetaCer,thetaTh[iPart],sigmaPID[iPart],kTRUE);
238     else                  height[iPart] = 0;
239     totalHeight +=height[iPart];
240     AliDebugClass(1,Form(" Particle %s with mass %f with height %f and thetaTH %f",AliPID::ParticleName(iPart),mass,height[iPart],thetaTh[iPart]));
241     AliDebugClass(1,Form(" partial height %15.14f total height %15.14f",height[iPart],totalHeight));
242   }
243   if(totalHeight<1e-5) {for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)richPID[iPart]=1.0/AliPID::kSPECIES;return;}
244   for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) richPID[iPart] = height[iPart]/totalHeight;
245   Int_t iPartNear = TMath::LocMax(AliPID::kSPECIES,richPID);
246   if(TMath::Abs(thetaCer-thetaTh[iPartNear])>5.*sigmaPID[iPartNear]) for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)richPID[iPart]=1.0/AliPID::kSPECIES;
247   //last line is to check if the nearest thetacerenkov to the teorethical one is within 5 sigma, otherwise no response (equal prob to every particle
248
249 }//CalcProb
250 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
251 void AliRICHTracker::EsdPrint()
252 {
253 // Reads a set of ESD files and print out some information
254 // Arguments: probCut - cut on probability 
255 //   Returns: none
256     
257   TFile *pFile=TFile::Open("AliESDs.root","read");              if(!pFile) {Printf("ERROR: AliESDs.root does not exist!");return;} 
258   TTree *pTr=(TTree*)pFile->Get("esdTree");                     if(!pTr)   {Printf("ERROR: AliESDs.root, no ESD tree inside!");return;} 
259   AliESD *pEsd=new AliESD;  pTr->SetBranchAddress("ESD", &pEsd); 
260   
261   Int_t iNevt=pTr->GetEntries();  Printf("This ESD contains %i events",iNevt);
262   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//ESD events loop
263     pTr->GetEvent(iEvt);
264     Int_t iNtracks=pEsd->GetNumberOfTracks(); Printf("ESD contains %i tracks created in Bz=%.2f Tesla",iNtracks,pEsd->GetMagneticField()/10.);
265     for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
266       AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
267       Float_t dx,dy;         pTrack->GetRICHdxdy(dx,dy);
268       Float_t theta,phi;     pTrack->GetRICHthetaPhi(theta,phi);
269       TString comment;
270       if(pTrack->GetRICHsignal()>0)
271         comment="OK";
272       else if(pTrack->GetRICHsignal()==kMipQdcCut)
273         comment="no enough QDC";
274       else if(pTrack->GetRICHsignal()==kMipDistCut)
275         comment="nearest cluster is too far";
276       else if(pTrack->GetRICHsignal()==-1)
277         comment="no intersection";            
278       Printf("Track %2i Q=%4.1f P=%.3f GeV RICH: ChamClus %7i Track-Mip=(%7.2f,%7.2f)=%5.2f cm ThetaCer %7.1f rad %s",
279                    iTrk,pTrack->GetSign(),pTrack->GetP(),
280                    pTrack->GetRICHcluster(),dx,dy,TMath::Sqrt(dx*dx+dy*dy),pTrack->GetRICHsignal(), comment.Data());
281     }//ESD tracks loop
282   }//ESD events loop
283   delete pEsd;  pFile->Close();//close AliESDs.root
284 }
285 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
286 void AliRICHTracker::MatrixPrint(Double_t probCut)
287 {
288 // Reads a set of 3  ESD files from current directory and prints out the matrix of probabilities to pion kaon or proton completely blindly withou nay assumption on the contents of files.
289 // Normally it implies that those 3 ESDs contain only particles of the same sort namly pions, kaons and protons in that order.  
290 // Arguments: probCut - cut on probability 
291 //   Returns: none
292   for(Int_t iFile=0;iFile<3;iFile++){
293     TFile *pFile=TFile::Open(Form("Esd%1i.root",iFile+1),"read"); if(!pFile) {Printf("ERROR: Esd%1i.root does not exist!",iFile+1);return;} 
294     TTree *pTr=(TTree*)pFile->Get("esdTree");                     if(!pTr)   {Printf("ERROR: Esd%1i.root, no ESD tree inside!",iFile+1);return;} 
295     AliESD *pEsd=new AliESD;  pTr->SetBranchAddress("ESD", &pEsd); 
296     Int_t iProtCnt=0,iKaonCnt=0,iPionCnt=0,iUnreconCnt=0,iTrkCnt=0; //counters
297     
298     for(Int_t iEvt=0;iEvt<pTr->GetEntries();iEvt++){//ESD events loop
299       pTr->GetEvent(iEvt);
300       iTrkCnt+=pEsd->GetNumberOfTracks();
301       for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//ESD tracks loop
302         AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
303         Float_t dx,dy;         pTrack->GetRICHdxdy(dx,dy);
304         Float_t theta,phi;     pTrack->GetRICHthetaPhi(theta,phi);
305         Double_t prob[5];       pTrack->GetRICHpid(prob);
306         if(pTrack->GetRICHsignal()>0){
307           if(prob[4]>probCut)                         iProtCnt++; 
308           if(prob[3]>probCut)                         iKaonCnt++;
309           if((prob[0]+prob[1]+prob[2])>probCut)       iPionCnt++;
310         } else
311           iUnreconCnt++;       
312       }//ESD tracks loop
313       
314     }//ESD events loop
315     Printf("Bz=%5.2f Events=%i Total tracks=%i No recognized tracks=%i Pion=%i Kaon=%i Proton=%i ProbCut=%.2f",
316         0.1*pEsd->GetMagneticField(),pTr->GetEntries(),iTrkCnt,iUnreconCnt,iPionCnt,iKaonCnt,iProtCnt,probCut);
317     delete pEsd;  pFile->Close();//close AliESDs.root
318   }//files loop
319 }
320 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++