1 #include "AliRICHTracker.h" //class header
3 #include "AliRICHRecon.h"
6 #include <TTree.h> //EsdPrint()
7 #include <TFile.h> //EsdPrint()
8 #include "AliRICHHelix.h"
11 #include <TParticle.h>
14 #include <TNtupleD.h> //RecWithStack();
15 #include <AliTrackPointArray.h> //GetTrackPoint()
16 #include <AliAlignObj.h> //GetTrackPoint()
17 ClassImp(AliRICHTracker)
18 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
19 AliRICHTracker::AliRICHTracker():AliTracker()
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;
27 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28 Bool_t AliRICHTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
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());
45 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
46 Int_t AliRICHTracker::LoadClusters(TTree *pCluTree)
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;
53 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
54 Int_t AliRICHTracker::PropagateBack(AliESD *pESD)
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"));
63 AliPID pid; // needed to retrive all the PID info
65 for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
66 AliESDtrack *pTrack = pESD->GetTrack(iTrk);// get next reconstructed track
67 Double_t mom[3], pos[3];
68 pTrack->GetPxPyPz(mom); TVector3 mom3(mom[0],mom[1],mom[2]);
69 pTrack->GetXYZ(pos); TVector3 pos3(pos[0],pos[1],pos[2]);
70 AliRICHHelix helix(pos3,mom3,(Int_t)pTrack->GetSign(),-0.1*GetBz()); //construct helix out of track running parameters
71 //Printf(" magnetic field %f charged %f\n",GetBz(),pTrack->GetSign()); helix.Print("Track");
72 Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());
73 if(!iChamber) continue; //no intersection with chambers, ignore this track go after the next one
75 //find MIP cluster candidate (closest to track intersection point cluster with large enough QDC)
76 Double_t dR=9999, dX=9999, dY=9999; //distance between track-PC intersection point and current cluster
77 Double_t mipDr=9999,mipDx=9999,mipDy=9999,mipX=9999,mipY=9999; //nearest cluster parameters
78 Int_t iMipId=-1; //index of this nearest cluster
79 for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
80 AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
81 if(pClu->Q()<AliRICHParam::QthMIP()) continue; //to low QDC, go after another one
82 pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY); //get distance for current cluster
83 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
84 }//clusters loop for intersected chamber
86 pTrack->SetRICHthetaPhi(helix.Ploc().Theta(),helix.Ploc().Phi()); //store track impact angles with respect to RICH planes
87 pTrack->SetRICHdxdy(mipDx,mipDy); //distance between track-PC intersection and closest cluster with Qdc>100
88 pTrack->SetRICHmipXY(mipX,mipY); //position of that closest cluster with Qdc>100
90 if(iMipId==-1) {pTrack->SetRICHsignal(kMipQdcCut); continue;} //no cluster with enough QDC found
91 if(mipDr>AliRICHParam::DmatchMIP()) {pTrack->SetRICHsignal(kMipDistCut); continue;} //closest cluster with enough carge is still too far
93 pTrack->SetRICHcluster(iMipId+1000000*iChamber); //set mip cluster index
94 pTrack->SetRICHsignal(recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId));//search for mean Cerenkov angle for this track
95 pTrack->SetRICHnclusters(iMipId); //on return iMipId is number of photon clusters accepted in reconstruction
97 AliDebug(1,Form("Ch=%i PC Intersection=(%5.2f,%5.2f) cm MIP cluster dist=(%5.2f,%5.2f)=%5.2f cm ThetaCkov=%f",
98 iChamber,helix.PosPc().X(),helix.PosPc().Y(), mipDx,mipDy,mipDr, pTrack->GetRICHsignal()));
100 //here comes PID calculations
101 if(pTrack->GetRICHsignal()>0) {
102 AliDebug(1,Form("Start to assign the probabilities"));
103 Double_t sigmaPID[AliPID::kSPECIES]; Double_t richPID[AliPID::kSPECIES];
104 for (Int_t iPid=0;iPid<AliPID::kSPECIES;iPid++){//PID loop
105 sigmaPID[iPid] = 0; fErrPar[iPid] = 0;
106 for(Int_t iCkov=0;iCkov<pRich->Clus(iChamber)->GetEntries();iCkov++){//Ckov candidates loop ????????????? somehting fomr AliRICHRecon must be
107 recon.SetPhotonIndex(iCkov);
108 if(recon.GetPhotonFlag() == 2) {
109 Double_t thetaCkov=0.6; //??????????????????
110 Double_t phiCkov=0.6; //??????????????????
111 Double_t thetaTrk=recon.GetTrackTheta();
112 Double_t phiTrk=(recon.GetPhiPoint()-recon.GetTrackPhi());
113 Double_t beta =pTrack->GetP()/TMath::Sqrt((pTrack->GetP()*pTrack->GetP()+pid.ParticleMass(iPid)*pid.ParticleMass(iPid)));
114 Double_t sigma2 = AliRICHParam::SigmaSinglePhotonFormula(thetaCkov,phiCkov,thetaTrk,phiTrk,beta);
115 if(sigma2>0) sigmaPID[iPid] += 1/sigma2;
117 }//Ckov candidates loop
119 if (sigmaPID[iPid]>0)
120 sigmaPID[iPid] *= (Double_t)(iMipId-recon.GetPhotBKG())/(Double_t)(iMipId); // n total phots, m are background...the sigma are scaled..
121 if(sigmaPID[iPid]>0) sigmaPID[iPid] = 1/TMath::Sqrt(sigmaPID[iPid])*0.001; // sigma from parametrization are in mrad...
122 else sigmaPID[iPid] = 0;
123 fErrPar[iPid]=sigmaPID[iPid];
124 AliDebug(1,Form("sigma for %s is %f rad",AliPID::ParticleName(iPid),sigmaPID[iPid]));
126 CalcProb(pTrack->GetRICHsignal(),pTrack->GetP(),sigmaPID,richPID);
127 pTrack->SetRICHpid(richPID);
128 AliDebug(1,Form("PROBABILITIES ---> %f - %f - %f - %f - %f",richPID[0],richPID[1],richPID[2],richPID[3],richPID[4]));
129 }//if(pTrack->GetRICHsignal())
131 AliDebug(1,"Stop pattern recognition");
132 return 0; // error code: 0=no error;
134 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
135 void AliRICHTracker::RecWithStack(TNtupleD *hn)
137 // 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.
138 // Arguments: hn- output ntuple where to store all variables
140 AliDebug(1,"Start.");
141 AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
143 // pRich->GetLoader()->GetRunLoader()->LoadHeader();
144 if(!pRich->GetLoader()->GetRunLoader()->TreeK()) pRich->GetLoader()->GetRunLoader()->LoadKinematics();
145 AliStack *pStack = pRich->GetLoader()->GetRunLoader()->Stack();
146 if(!pStack) {AliDebug(1,Form("No STACK found in AliRoot"));return;}
147 Int_t iNtracks=pStack->GetNtrack();
148 AliDebug(1,Form(" Start reconstruction with %i track(s) from Stack",iNtracks));
152 Double_t b=GetFieldMap()->SolenoidField()/10;// magnetic field in Tesla
153 AliDebug(1,Form("Start with simulated %i tracks in %f Tesla field",iNtracks,b));
154 TVector3 x0(0,0,0); TVector3 p0(0,0,0);//tmp storage for AliRICHHelix
157 if(pRich->GetLoader()->LoadRecPoints()) {AliDebug(1,Form("No clusters found in RICH"));return;}
158 pRich->GetLoader()->TreeR()->GetEntry(0);
161 for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//stack particles loop
162 TParticle *pParticle = pStack->Particle(iTrackN);
163 if(!pParticle) {AliDebug(1,Form("Not a valid TParticle pointer. Track skipped"));continue;}
164 AliDebug(1,Form(" PDG code : %i",pParticle->GetPdgCode()));
166 // problem of PDG code of some extra particles to be solved!!!!!!!!!
168 // found problem! Look in TRD directory : codes from Fluka are :
170 // if ((pdg_code == 10010020) ||
171 // (pdg_code == 10010030) ||
172 // (pdg_code == 50000050) ||
173 // (pdg_code == 50000051) ||
174 // (pdg_code == 10020040)) {
176 if(pParticle->GetPdgCode()>=50000050||pParticle->GetPdgCode()==0||pParticle->GetPdgCode()>10000) {AliDebug(1,Form("A photon as track... Track skipped"));continue;}
177 if(!pParticle->GetPDG()) continue;
179 // to be updated for us!!
181 AliDebug(1,Form("Track %i is a %s with charge %i and momentum %f",
182 iTrackN,pParticle->GetPDG()->GetName(),Int_t(pParticle->GetPDG()->Charge()),pParticle->P()));
183 // if(pParticle->GetMother(0)!=-1) continue; //consider only primaries
184 if(pParticle->GetPDG()->Charge()==0||TMath::Abs(Int_t(pParticle->GetPDG()->Charge()))!=3) continue; //to avoid photons from stack...
185 hnvec[0]=pParticle->P();
186 hnvec[1]=pParticle->GetPDG()->Charge();
188 p0.SetMagThetaPhi(pParticle->P(),pParticle->Theta(),pParticle->Phi()); x0.SetXYZ(pParticle->Vx(),pParticle->Vy(),pParticle->Vz());
189 AliRICHHelix helix(x0,p0,TMath::Sign(1,(Int_t)pParticle->GetPDG()->Charge()),b);
190 Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());
191 if(!iChamber) continue;// no intersection with RICH found
193 hnvec[2]=helix.Ploc().Theta();
194 hnvec[3]=helix.Ploc().Phi();
195 hnvec[4]=helix.PosPc().X();
196 hnvec[5]=helix.PosPc().Y();
198 Double_t dX,dY,dR,dRmip=9999; //min distance between clusters and track position on PC
199 Int_t iMipId=-1; //index of that min distance cluster
200 for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
201 AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
202 pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY);//ditance between current cluster and helix intersection with PC
203 if(dR<dRmip){dRmip=dR; hnvec[6]=pClu->X();hnvec[7]=pClu->Y();hnvec[8]=pClu->Q();
204 iMipId=1000000*iChamber+iClu;}//find cluster nearest to the track
206 }//clusters loop for intersected chamber
209 hnvec[9]=recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId); //search for mean Cerenkov angle for this track
210 hnvec[10]=iMipId;//on return from ThetaCerenkov() contains number of photon candidates accepted
211 hnvec[11]=(Double_t)iMipId;
212 hnvec[12]=(Double_t)iChamber;
213 hnvec[13]=(Double_t)pParticle->GetPdgCode();
214 if(hn) hn->Fill(hnvec);
215 AliDebug(1,Form("FINAL Theta Cerenkov=%f",hnvec[9]));
216 }//stack particles loop
218 pRich->GetLoader()->UnloadRecPoints();
221 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
222 void AliRICHTracker::CalcProb(Double_t thetaCer,Double_t pmod, Double_t *sigmaPID, Double_t *richPID)
224 // Calculates probability to be a electron-muon-pion-kaon-proton
225 // from the given Cerenkov angle and momentum assuming no initial particle composition
226 // (i.e. apriory probability to be the particle of the given sort is the same for all sorts)
227 Double_t height[AliPID::kSPECIES];Double_t totalHeight=0;
228 Double_t thetaTh[AliPID::kSPECIES];
229 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
231 Double_t mass = AliRICHParam::fgMass[iPart];
232 Double_t refIndex=AliRICHParam::Instance()->IdxC6F14(AliRICHParam::EckovMean());
233 Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(refIndex*pmod);
235 if(cosThetaTh>=1) continue;
236 thetaTh[iPart] = TMath::ACos(cosThetaTh);
237 // Double_t sinThetaThNorm = TMath::Sin(thetaTh)/TMath::Sqrt(1-1/(refIndex*refIndex));
238 // Double_t sigmaThetaTh = (0.014*(1/sinThetaThNorm-1) + 0.0043)*1.25;
239 // height[iPart] = TMath::Gaus(thetaCer,thetaTh,sigmaThetaTh);
240 if(sigmaPID[iPart]>0) height[iPart] = TMath::Gaus(thetaCer,thetaTh[iPart],sigmaPID[iPart],kTRUE);
241 else height[iPart] = 0;
242 totalHeight +=height[iPart];
243 AliDebugClass(1,Form(" Particle %s with mass %f with height %f and thetaTH %f",AliPID::ParticleName(iPart),mass,height[iPart],thetaTh[iPart]));
244 AliDebugClass(1,Form(" partial height %15.14f total height %15.14f",height[iPart],totalHeight));
246 if(totalHeight<1e-5) {for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)richPID[iPart]=1.0/AliPID::kSPECIES;return;}
247 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) richPID[iPart] = height[iPart]/totalHeight;
248 Int_t iPartNear = TMath::LocMax(AliPID::kSPECIES,richPID);
249 if(TMath::Abs(thetaCer-thetaTh[iPartNear])>5.*sigmaPID[iPartNear]) for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)richPID[iPart]=1.0/AliPID::kSPECIES;
250 //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
253 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
254 void AliRICHTracker::EsdPrint()
256 // Reads a set of ESD files and print out some information
257 // Arguments: probCut - cut on probability
260 TFile *pFile=TFile::Open("AliESDs.root","read"); if(!pFile) {Printf("ERROR: AliESDs.root does not exist!");return;}
261 TTree *pTr=(TTree*)pFile->Get("esdTree"); if(!pTr) {Printf("ERROR: AliESDs.root, no ESD tree inside!");return;}
262 AliESD *pEsd=new AliESD; pTr->SetBranchAddress("ESD", &pEsd);
264 Int_t iNevt=pTr->GetEntries(); Printf("This ESD contains %i events",iNevt);
265 for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//ESD events loop
267 Int_t iNtracks=pEsd->GetNumberOfTracks(); Printf("ESD contains %i tracks created in Bz=%.2f Tesla",iNtracks,pEsd->GetMagneticField()/10.);
268 for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
269 AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
270 Float_t dx,dy; pTrack->GetRICHdxdy(dx,dy);
271 Float_t theta,phi; pTrack->GetRICHthetaPhi(theta,phi);
273 if(pTrack->GetRICHsignal()>0)
275 else if(pTrack->GetRICHsignal()==kMipQdcCut)
276 comment="no enough QDC";
277 else if(pTrack->GetRICHsignal()==kMipDistCut)
278 comment="nearest cluster is too far";
279 else if(pTrack->GetRICHsignal()==-1)
280 comment="no intersection";
281 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",
282 iTrk,pTrack->GetSign(),pTrack->GetP(),
283 pTrack->GetRICHcluster(),dx,dy,TMath::Sqrt(dx*dx+dy*dy),pTrack->GetRICHsignal(), comment.Data());
286 delete pEsd; pFile->Close();//close AliESDs.root
288 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
289 void AliRICHTracker::MatrixPrint(Double_t probCut)
291 // 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.
292 // Normally it implies that those 3 ESDs contain only particles of the same sort namly pions, kaons and protons in that order.
293 // Arguments: probCut - cut on probability
295 for(Int_t iFile=0;iFile<3;iFile++){
296 TFile *pFile=TFile::Open(Form("Esd%1i.root",iFile+1),"read"); if(!pFile) {Printf("ERROR: Esd%1i.root does not exist!",iFile+1);return;}
297 TTree *pTr=(TTree*)pFile->Get("esdTree"); if(!pTr) {Printf("ERROR: Esd%1i.root, no ESD tree inside!",iFile+1);return;}
298 AliESD *pEsd=new AliESD; pTr->SetBranchAddress("ESD", &pEsd);
299 Int_t iProtCnt=0,iKaonCnt=0,iPionCnt=0,iUnreconCnt=0,iTrkCnt=0; //counters
301 for(Int_t iEvt=0;iEvt<pTr->GetEntries();iEvt++){//ESD events loop
303 iTrkCnt+=pEsd->GetNumberOfTracks();
304 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//ESD tracks loop
305 AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
306 Float_t dx,dy; pTrack->GetRICHdxdy(dx,dy);
307 Float_t theta,phi; pTrack->GetRICHthetaPhi(theta,phi);
308 Double_t prob[5]; pTrack->GetRICHpid(prob);
309 if(pTrack->GetRICHsignal()>0){
310 if(prob[4]>probCut) iProtCnt++;
311 if(prob[3]>probCut) iKaonCnt++;
312 if((prob[0]+prob[1]+prob[2])>probCut) iPionCnt++;
318 Printf("Bz=%5.2f Events=%i Total tracks=%i No recognized tracks=%i Pion=%i Kaon=%i Proton=%i ProbCut=%.2f",
319 0.1*pEsd->GetMagneticField(),pTr->GetEntries(),iTrkCnt,iUnreconCnt,iPionCnt,iKaonCnt,iProtCnt,probCut);
320 delete pEsd; pFile->Close();//close AliESDs.root
323 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++