1a86e3b4037d081ab53fe3ea10a8fa0aac1435d6
[u/mrichter/AliRoot.git] / RICH / AliRICHTracker.cxx
1 #include "AliRICHTracker.h" //class header
2 #include "AliRICH.h"
3 #include "AliRICHRecon.h"
4 #include "AliRICHParam.h"
5 #include <AliESD.h>
6 #include <TGeoManager.h>     //EsdQA()
7 #include <TVector3.h>
8 #include <TTree.h>          //EsdQA() 
9 #include <TFile.h>          //EsdQA()  
10 #include <TProfile.h>   //EsdQA() 
11 #include "AliRICHHelix.h"
12 #include <AliMagF.h>
13 #include <AliStack.h>
14 #include <TParticle.h>
15 #include <TMath.h>
16 #include <AliRun.h>
17 #include <TNtupleD.h>            //RecWithStack();
18 #include <AliTrackPointArray.h>  //GetTrackPoint()
19 #include <AliAlignObj.h>         //GetTrackPoint()
20 #include <TH1F.h>                //EsdQA()  
21 #include <TH2F.h>                //EsdQA()  
22 #include <TCanvas.h>             //EsdQA()  
23 ClassImp(AliRICHTracker)
24 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
25 AliRICHTracker::AliRICHTracker():AliTracker()
26 {
27 // AliRICHTracker is created from AliReconstraction::Run() which invokes AliReconstraction::CreateTrackers() 
28 // which in turn invokes AliRICHReconstructor::CreateTracker(). 
29 // Note that this is done just once per session before AliReconstruction::Run() goes to events loop.
30   AliRICHParam::Instance()->CdbRead(0,0); 
31   for(Int_t i=0;i<5;i++)fErrPar[i]=0;
32 }
33 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
34 Bool_t AliRICHTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
35 {
36 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
37 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
38 // Arguments: idx- cluster index which is stored by RICH in AliESDtrack
39 //            point- reference to the object where to store the point     
40 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
41   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
42   Int_t iCham=idx/1000000;
43   Int_t iClu=idx%1000000;
44   point.SetVolumeID(AliAlignObj::LayerToVolUID(AliAlignObj::kRICH,iCham-1));//layer and chamber number
45   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));  
46   AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iCham)->UncheckedAt(iClu);//get pointer to cluster
47   TVector3 mars=AliRICHParam::Instance()->Lors2Mars(iCham,pClu->X(),pClu->Y());
48   point.SetXYZ(mars.X(),mars.Y(),mars.Z());
49   return kTRUE;
50 }
51 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
52 Int_t AliRICHTracker::LoadClusters(TTree *pCluTree)
53 {
54 // Interface callback methode invoked from AliReconstruction::RunTracking() to load RICH clusters for RICH
55 // Arguments: pCluTree- pointer to clusters tree got by AliRICHLoader::LoadRecPoints("read") then AliRICHLoader::TreeR()
56 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
57   AliDebug(1,"Start.");  pCluTree->GetEntry(0);  AliDebug(1,"Stop."); return 0;
58 }
59 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 Int_t AliRICHTracker::PropagateBack(AliESD *pESD)
61 {
62 // Interface callback methode invoked by AliRecontruction::RunTracking() during tracking after TOF. It's done just once per event
63 // Arguments: pESD - pointer to Event Summary Data class instance which contains a list of tracks
64 //   Returns: error code, 0 if no errors   
65   Int_t iNtracks=pESD->GetNumberOfTracks();
66   AliDebug(1,Form("Start with %i tracks",iNtracks));
67   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));  
68   AliRICHRecon recon;
69   AliPID pid; // needed to retrive all the PID info
70    
71   for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
72     AliESDtrack *pTrack = pESD->GetTrack(iTrk);// get next reconstructed track    
73     Double_t mom[3], pos[3];
74     pTrack->GetPxPyPz(mom); TVector3 mom3(mom[0],mom[1],mom[2]);
75     pTrack->GetXYZ(pos); TVector3 pos3(pos[0],pos[1],pos[2]);
76     AliRICHHelix helix(pos3,mom3,(Int_t)pTrack->GetSign(),-0.1*GetBz()); //construct helix out of track running parameters  
77      //Printf(" magnetic field %f charged %f\n",GetBz(),pTrack->GetSign()); helix.Print("Track");
78     Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());    
79     if(!iChamber) continue;                                         //no intersection with chambers, ignore this track go after the next one
80   
81     //find MIP cluster candidate (closest to track intersection point cluster with large enough QDC)
82     Double_t    dR=9999,   dX=9999,   dY=9999;                     //distance between track-PC intersection point and current cluster
83     Double_t mipDr=9999,mipDx=9999,mipDy=9999,mipX=9999,mipY=9999; //nearest cluster parameters
84     Int_t   iMipId=-1;                         //index of this nearest cluster
85     for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
86       AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
87       if(pClu->Q()<AliRICHParam::QthMIP()) continue; //to low QDC, go after another one
88       pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY); //get distance for current cluster
89       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
90     }//clusters loop for intersected chamber
91
92     pTrack->SetRICHthetaPhi(helix.Ploc().Theta(),helix.Ploc().Phi()); //store track impact angles with respect to RICH planes
93     pTrack->SetRICHdxdy(mipDx,mipDy);                                 //distance between track-PC intersection and closest cluster with Qdc>100
94     pTrack->SetRICHmipXY(mipX,mipY);                                  //position of that closest cluster with Qdc>100
95     
96     if(iMipId==-1)                        {pTrack->SetRICHsignal(kMipQdcCut);  continue;} //no cluster with enough QDC found
97     if(mipDr>AliRICHParam::DmatchMIP())   {pTrack->SetRICHsignal(kMipDistCut); continue;} //closest cluster with enough carge is still too far 
98   
99     pTrack->SetRICHcluster(iMipId+1000000*iChamber);                                //set mip cluster index
100     Int_t nphots =0;
101     pTrack->SetRICHsignal(recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),nphots));//search for mean Cerenkov angle for this track
102     pTrack->SetRICHchi2(recon.GetRingSigma2());
103     pTrack->SetRICHnclusters(nphots);                                               //on return nphots is number of photon clusters accepted in reconstruction
104
105     AliDebug(1,Form("Ch=%i PC Intersection=(%5.2f,%5.2f) cm MIP cluster dist=(%5.2f,%5.2f)=%5.2f cm ThetaCkov=%f",
106                      iChamber,helix.PosPc().X(),helix.PosPc().Y(),            mipDx,mipDy,mipDr,     pTrack->GetRICHsignal()));
107     
108 //here comes PID calculations    
109 //    CalcProb(pTrack->GetRICHsignal(),pTrack->GetP(),sigmaPID,richPID);
110   }//ESD tracks loop
111   AliDebug(1,"Stop pattern recognition");
112   return 0; // error code: 0=no error;
113 }//PropagateBack()
114 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
115 void AliRICHTracker::RecWithStack(TNtupleD *hn)
116 {
117 // 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.
118 // Arguments: hn- output ntuple where to store all variables
119 //   Returns: none       
120   AliDebug(1,"Start.");  
121   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
122   
123 //  pRich->GetLoader()->GetRunLoader()->LoadHeader();
124   if(!pRich->GetLoader()->GetRunLoader()->TreeK()) pRich->GetLoader()->GetRunLoader()->LoadKinematics();
125   AliStack *pStack =   pRich->GetLoader()->GetRunLoader()->Stack();
126   if(!pStack) {AliDebug(1,Form("No STACK found in AliRoot"));return;}
127   Int_t iNtracks=pStack->GetNtrack();
128   AliDebug(1,Form(" Start reconstruction with %i track(s) from Stack",iNtracks));
129   
130   Double_t hnvec[20];
131   
132   Double_t b=GetFieldMap()->SolenoidField()/10;// magnetic field in Tesla
133   AliDebug(1,Form("Start with simulated %i tracks in %f Tesla field",iNtracks,b));
134   TVector3 x0(0,0,0); TVector3 p0(0,0,0);//tmp storage for AliRICHHelix
135   
136
137   if(pRich->GetLoader()->LoadRecPoints()) {AliDebug(1,Form("No clusters found in RICH"));return;}
138   pRich->GetLoader()->TreeR()->GetEntry(0);
139
140   AliRICHRecon recon;
141   for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//stack particles loop
142     TParticle *pParticle = pStack->Particle(iTrackN);
143     if(!pParticle) {AliDebug(1,Form("Not a valid TParticle pointer. Track skipped"));continue;}
144     AliDebug(1,Form(" PDG code : %i",pParticle->GetPdgCode()));
145 //
146 // problem of PDG code of some extra particles to be solved!!!!!!!!!
147 //
148 // found problem! Look in TRD directory : codes from Fluka are :
149 //
150 //    if ((pdg_code == 10010020) ||
151 //        (pdg_code == 10010030) ||
152 //        (pdg_code == 50000050) ||
153 //        (pdg_code == 50000051) ||
154 //        (pdg_code == 10020040)) {
155 //
156     if(pParticle->GetPdgCode()>=50000050||pParticle->GetPdgCode()==0||pParticle->GetPdgCode()>10000) {AliDebug(1,Form("A photon as track... Track skipped"));continue;}
157     if(!pParticle->GetPDG()) continue;
158 //
159 // to be updated for us!!
160 //
161     AliDebug(1,Form("Track %i is a %s with charge %i and momentum %f",
162             iTrackN,pParticle->GetPDG()->GetName(),Int_t(pParticle->GetPDG()->Charge()),pParticle->P()));
163 //    if(pParticle->GetMother(0)!=-1) continue; //consider only primaries
164     if(pParticle->GetPDG()->Charge()==0||TMath::Abs(Int_t(pParticle->GetPDG()->Charge()))!=3) continue; //to avoid photons from stack...
165     hnvec[0]=pParticle->P();
166     hnvec[1]=pParticle->GetPDG()->Charge();
167     
168     p0.SetMagThetaPhi(pParticle->P(),pParticle->Theta(),pParticle->Phi());   x0.SetXYZ(pParticle->Vx(),pParticle->Vy(),pParticle->Vz());
169     AliRICHHelix helix(x0,p0,TMath::Sign(1,(Int_t)pParticle->GetPDG()->Charge()),b);   
170     Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());        
171     if(!iChamber) continue;// no intersection with RICH found
172     
173     hnvec[2]=helix.Ploc().Theta();
174     hnvec[3]=helix.Ploc().Phi();
175     hnvec[4]=helix.PosPc().X();
176     hnvec[5]=helix.PosPc().Y();
177     
178     Double_t dX,dY,dR,dRmip=9999;      //min distance between clusters and track position on PC 
179     Int_t iMipId=-1;       //index of that min distance cluster 
180     for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
181       AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
182       pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY);//ditance between current cluster and helix intersection with PC
183       if(dR<dRmip){dRmip=dR; hnvec[6]=pClu->X();hnvec[7]=pClu->Y();hnvec[8]=pClu->Q();
184                                                   iMipId=1000000*iChamber+iClu;}//find cluster nearest to the track 
185       
186     }//clusters loop for intersected chamber
187     
188     
189     hnvec[9]=recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId); //search for mean Cerenkov angle for this track
190     hnvec[10]=iMipId;//on return from ThetaCerenkov() contains number of photon candidates accepted
191     hnvec[11]=(Double_t)iMipId;
192     hnvec[12]=(Double_t)iChamber;
193     hnvec[13]=(Double_t)pParticle->GetPdgCode();
194     if(hn) hn->Fill(hnvec);
195     AliDebug(1,Form("FINAL Theta Cerenkov=%f",hnvec[9]));
196   }//stack particles loop
197   
198   pRich->GetLoader()->UnloadRecPoints();
199   AliDebug(1,"Stop.");  
200 }//RecWithStack
201 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
202 void AliRICHTracker::EsdQA(Bool_t isPrint)
203 {
204 // Reads a set of ESD files and print out or plot some information
205 // Arguments: isPrint is a flag to choose between printing (isPrint = kTRUE) and plotting (isPrint = kFALSE)
206 //   Returns: none
207     
208   TFile *pFile=TFile::Open("AliESDs.root","read");              if(!pFile) {Printf("ERROR: AliESDs.root does not exist!");return;} 
209   TTree *pTr=(TTree*)pFile->Get("esdTree");                     if(!pTr)   {Printf("ERROR: AliESDs.root, no ESD tree inside!");return;} 
210   AliESD *pEsd=new AliESD;  pTr->SetBranchAddress("ESD", &pEsd); 
211   
212   Int_t iNevt=pTr->GetEntries();  Printf("This ESD contains %i events",iNevt);
213    //AliRICHParam *pPar;
214   
215   TH1D *pDx=0,*pDy=0,*pProbEl=0,*pProbMu=0,*pProbPi=0,*pProbKa=0,*pProbPr=0; TH2F *pThP=0; TProfile *pChiTh=0;
216   if(!isPrint){
217     TH1::AddDirectory(kFALSE);    
218     pDx    =new TH1D("dX"  ,"distance between Xmip and Xtrack;cm",300,-1.5,1.5);
219     pDy    =new TH1D("dY"  ,"distance between Ymip and Ytrack;cm",300,-1.5,1.5);
220     pThP   =new TH2F("tvsp","#theta_{Ckov} radian;P GeV"              ,65 ,-0.5,6.0,75,0,0.75); pThP->SetStats(0);
221     pProbPi=new TH1D("RichProbPion","HMPID PID probability for  e #mu #pi"   ,101,0.05,1.05); pProbPi->SetLineColor(kRed);
222     pProbEl=new TH1D("RichProbEle" ,""   ,101,0.05,1.05);                                     pProbEl->SetLineColor(kGreen); 
223     pProbMu=new TH1D("RichProbMuon" ,""  ,101,0.05,1.05);                                     pProbMu->SetLineColor(kBlue);
224     pProbKa=new TH1D("RichProbKaon","HMPID PID probability for K"     ,101,0.05,1.05);
225     pProbPr=new TH1D("RichProbProton","HMPID PID probability for p"   ,101,0.05,1.05);
226     pChiTh =new TProfile("RichChiTh","#chi^{2};#theta_{C}"            ,80 ,0,0.8 , -2,2);
227     
228     //if(!gGeoManager)TGeoManager::Import("geometry.root");
229     //pPar = AliRICHParam::Instance();
230   }
231   
232   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//ESD events loop
233     pTr->GetEvent(iEvt);
234     Int_t iNtracks=pEsd->GetNumberOfTracks(); Printf("ESD contains %i tracks created in Bz=%.2f Tesla",iNtracks,pEsd->GetMagneticField()/10.);
235     for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
236       AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
237       Float_t dx,dy;         pTrack->GetRICHdxdy(dx,dy);
238       Float_t theta,phi;     pTrack->GetRICHthetaPhi(theta,phi);
239       TString comment;
240       if(isPrint){
241         if(pTrack->GetRICHsignal()>0)                         comment="OK";
242         else if(pTrack->GetRICHsignal()==kMipQdcCut)          comment="no enough QDC";
243         else if(pTrack->GetRICHsignal()==kMipDistCut)         comment="nearest cluster is too far";
244         else if(pTrack->GetRICHsignal()==-1)                  comment="no intersection";
245         Double_t pid[5]; pTrack->GetRICHpid(pid);           
246         Printf("Tr %2i Q=%4.1f P=%.3f GeV Tr-Mip=%5.2f cm Th=%7.3f rad Prob : e=%.4f mu=%.4f pi=%.4f K=%.4f p=%.4f %s" ,
247                      iTrk,pTrack->GetSign(),pTrack->GetP(),TMath::Sqrt(dx*dx+dy*dy),pTrack->GetRICHsignal(), 
248                      pid[0],pid[1],pid[2],pid[3],pid[4], comment.Data());
249       }else{//collect hists
250     
251         pDx->Fill(dx);  pDy->Fill(dy);
252         pThP->Fill(pTrack->GetP(),pTrack->GetRICHsignal());
253         pChiTh->Fill(pTrack->GetRICHsignal(),pTrack->GetRICHchi2());
254         Double_t prob[5]; pTrack->GetRICHpid(prob);
255         pProbEl->Fill(prob[0]); pProbMu->Fill(prob[1]); pProbPi->Fill(prob[2]);  pProbKa->Fill(prob[3]);  pProbPr->Fill(prob[4]); 
256       }//if plot
257     }//ESD tracks loop
258   }//ESD events loop
259   delete pEsd;  pFile->Close();//close AliESDs.root
260   if(!isPrint){
261     TCanvas *pC=new TCanvas("c","Quality",1200,1500); pC->Divide(2,2);
262     TF1 *pPion = new TF1("RICHtheor","acos(sqrt(x*x+[0]*[0])/(x*[1]))",1.2,6);  pPion->SetLineWidth(1);pPion->SetParameter(1,1.292);
263     AliPID ppp;                   pPion->SetParameter(0,AliPID::ParticleMass(AliPID::kPion))  ; pPion->SetLineColor(kRed);
264     TF1 *pKaon = (TF1*)pPion->Clone();  pKaon->SetParameter(0,AliPID::ParticleMass(AliPID::kKaon))  ; pKaon->SetLineColor(kGreen);
265     TF1 *pProt = (TF1*)pPion->Clone();  pProt->SetParameter(0,AliPID::ParticleMass(AliPID::kProton)); pProt->SetLineColor(kBlue);
266     
267     
268     pC->cd(1);        pDx->Draw();
269     pC->cd(2);        pDy->Draw();
270     pC->cd(3);        pThP->Draw();       pPion->Draw("same"); pKaon->Draw("same"); pProt->Draw("same");
271     pC->cd(4);        pChiTh->Draw();
272     
273     TCanvas *pC2=new TCanvas("c2","Quality 2",1200,1500); pC2->Divide(2,2);
274     pC2->cd(1);        pProbPi->Draw(); pProbMu->Draw("same"); pProbEl->Draw("same");
275     pC2->cd(2);        pProbKa->Draw();
276     pC2->cd(3);        pProbPr->Draw();
277     
278     
279   }
280 }
281 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
282 void AliRICHTracker::MatrixPrint(Double_t probCut)
283 {
284 // 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.
285 // Normally it implies that those 3 ESDs contain only particles of the same sort namly pions, kaons and protons in that order.  
286 // Arguments: probCut - cut on probability 
287 //   Returns: none
288   for(Int_t iFile=0;iFile<3;iFile++){
289     TFile *pFile=TFile::Open(Form("Esd%1i.root",iFile+1),"read"); if(!pFile) {Printf("ERROR: Esd%1i.root does not exist!",iFile+1);return;} 
290     TTree *pTr=(TTree*)pFile->Get("esdTree");                     if(!pTr)   {Printf("ERROR: Esd%1i.root, no ESD tree inside!",iFile+1);return;} 
291     AliESD *pEsd=new AliESD;  pTr->SetBranchAddress("ESD", &pEsd); 
292     Int_t iProtCnt=0,iKaonCnt=0,iPionCnt=0,iUnreconCnt=0,iTrkCnt=0; //counters
293     
294     for(Int_t iEvt=0;iEvt<pTr->GetEntries();iEvt++){//ESD events loop
295       pTr->GetEvent(iEvt);
296       iTrkCnt+=pEsd->GetNumberOfTracks();
297       for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//ESD tracks loop
298         AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
299         Float_t dx,dy;         pTrack->GetRICHdxdy(dx,dy);
300         Float_t theta,phi;     pTrack->GetRICHthetaPhi(theta,phi);
301         Double_t prob[5];       pTrack->GetRICHpid(prob);
302         if(pTrack->GetRICHsignal()>0){
303           if(prob[4]>probCut)                         iProtCnt++; 
304           if(prob[3]>probCut)                         iKaonCnt++;
305           if((prob[0]+prob[1]+prob[2])>probCut)       iPionCnt++;
306         } else
307           iUnreconCnt++;       
308       }//ESD tracks loop
309       
310     }//ESD events loop
311     Printf("Bz=%5.2f Events=%i Total tracks=%i No recognized tracks=%i Pion=%i Kaon=%i Proton=%i ProbCut=%.2f",
312         0.1*pEsd->GetMagneticField(),pTr->GetEntries(),iTrkCnt,iUnreconCnt,iPionCnt,iKaonCnt,iProtCnt,probCut);
313     delete pEsd;  pFile->Close();//close AliESDs.root
314   }//files loop
315 }
316 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++