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"));
64 for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
65 AliESDtrack *pTrack = pESD->GetTrack(iTrk);// get next reconstructed track
66 AliRICHHelix helix(pTrack->X3(),pTrack->P3(),(Int_t)pTrack->GetSign(),-GetBz()/10); //construct helix out of track running parameters
67 //Printf(" magnetic field %f charged %f\n",GetBz(),pTrack->GetSign()); helix.Print("Track");
68 Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());
69 if(!iChamber) continue; //no intersection with chambers, ignore this track go after the next one
71 //find MIP cluster candidate (closest to track intersection point cluster with large enough QDC)
72 Double_t dR=9999, dX=9999, dY=9999; //distance between track-PC intersection point and current cluster
73 Double_t dRmip=9999,dXmip=9999,dYmip=9999; //distance between track-PC intersection point and nearest cluster
74 Int_t iMipId=-1; //index of this nearest cluster
75 for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
76 AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
77 if(pClu->Q()<AliRICHParam::QthMIP()) continue; //to low QDC, go after another one
78 pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY); //get distance for current cluster
79 if(dR<dRmip){iMipId=iClu; dRmip=dR; dXmip=dX; dYmip=dY; } //current cluster is closer, overwrite data for min cluster
80 }//clusters loop for intersected chamber
82 pTrack->SetRICHthetaPhi(helix.Ploc().Theta(),helix.Ploc().Phi()); //store track impact angles with respect to RICH planes
83 pTrack->SetRICHdxdy(dXmip,dYmip); //distance between track-PC intersection and closest cluster with Qdc>100
85 if(iMipId==-1) {pTrack->SetRICHsignal(kMipQdcCut); continue;} //no cluster with enough QDC found
86 if(dRmip>AliRICHParam::DmatchMIP()) {pTrack->SetRICHsignal(kMipDistCut); continue;} //closest cluster with enough carge is still too far
88 pTrack->SetRICHcluster(iMipId+1000000*iChamber); //set mip cluster index
89 pTrack->SetRICHsignal(recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId));//search for mean Cerenkov angle for this track
90 pTrack->SetRICHnclusters(iMipId); //on return iMipId is number of photon clusters accepted in reconstruction
92 AliDebug(1,Form("Ch=%i PC Intersection=(%5.2f,%5.2f) cm MIP cluster dist=(%5.2f,%5.2f)=%5.2f cm ThetaCkov=%f",
93 iChamber,helix.PosPc().X(),helix.PosPc().Y(), dXmip,dYmip,dRmip, pTrack->GetRICHsignal()));
95 //here comes PID calculations
96 if(pTrack->GetRICHsignal()>0) {
97 AliDebug(1,Form("Start to assign the probabilities"));
98 Double_t sigmaPID[AliPID::kSPECIES];
99 Double_t richPID[AliPID::kSPECIES];
100 for (Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) {
103 for(Int_t iphot=0;iphot<pRich->Clus(iChamber)->GetEntries();iphot++) {
104 recon.SetPhotonIndex(iphot);
105 if(recon.GetPhotonFlag() == 2) {
106 Double_t theta_g=recon.GetTrackTheta();
107 Double_t phi_g=(recon.GetPhiPoint()-recon.GetTrackPhi());
108 Double_t sigma2 = AliRICHParam::Instance()->SigmaSinglePhoton(iPart,pTrack->GetP(),theta_g,phi_g).Mag2();
109 if(sigma2>0) sigmaPID[iPart] += 1/sigma2;
112 if (sigmaPID[iPart]>0)
113 sigmaPID[iPart] *= (Double_t)(iMipId-recon.GetPhotBKG())/(Double_t)(iMipId); // n total phots, m are background...the sigma are scaled..
114 if(sigmaPID[iPart]>0) sigmaPID[iPart] = 1/TMath::Sqrt(sigmaPID[iPart])*0.001; // sigma from parametrization are in mrad...
115 else sigmaPID[iPart] = 0;
116 fErrPar[iPart]=sigmaPID[iPart];
117 AliDebug(1,Form("sigma for %s is %f rad",AliPID::ParticleName(iPart),sigmaPID[iPart]));
119 CalcProb(pTrack->GetRICHsignal(),pTrack->GetP(),sigmaPID,richPID);
120 pTrack->SetRICHpid(richPID);
121 AliDebug(1,Form("PROBABILITIES ---> %f - %f - %f - %f - %f",richPID[0],richPID[1],richPID[2],richPID[3],richPID[4]));
124 AliDebug(1,"Stop pattern recognition");
125 return 0; // error code: 0=no error;
127 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
128 void AliRICHTracker::RecWithStack(TNtupleD *hn)
130 // 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.
131 // Arguments: hn- output ntuple where to store all variables
133 AliDebug(1,"Start.");
134 AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
136 // pRich->GetLoader()->GetRunLoader()->LoadHeader();
137 if(!pRich->GetLoader()->GetRunLoader()->TreeK()) pRich->GetLoader()->GetRunLoader()->LoadKinematics();
138 AliStack *pStack = pRich->GetLoader()->GetRunLoader()->Stack();
139 if(!pStack) {AliDebug(1,Form("No STACK found in AliRoot"));return;}
140 Int_t iNtracks=pStack->GetNtrack();
141 AliDebug(1,Form(" Start reconstruction with %i track(s) from Stack",iNtracks));
145 Double_t b=GetFieldMap()->SolenoidField()/10;// magnetic field in Tesla
146 AliDebug(1,Form("Start with simulated %i tracks in %f Tesla field",iNtracks,b));
147 TVector3 x0(0,0,0); TVector3 p0(0,0,0);//tmp storage for AliRICHHelix
150 if(pRich->GetLoader()->LoadRecPoints()) {AliDebug(1,Form("No clusters found in RICH"));return;}
151 pRich->GetLoader()->TreeR()->GetEntry(0);
154 for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//stack particles loop
155 TParticle *pParticle = pStack->Particle(iTrackN);
156 if(!pParticle) {AliDebug(1,Form("Not a valid TParticle pointer. Track skipped"));continue;}
157 AliDebug(1,Form(" PDG code : %i",pParticle->GetPdgCode()));
159 // problem of PDG code of some extra particles to be solved!!!!!!!!!
161 // found problem! Look in TRD directory : codes from Fluka are :
163 // if ((pdg_code == 10010020) ||
164 // (pdg_code == 10010030) ||
165 // (pdg_code == 50000050) ||
166 // (pdg_code == 50000051) ||
167 // (pdg_code == 10020040)) {
169 if(pParticle->GetPdgCode()>=50000050||pParticle->GetPdgCode()==0||pParticle->GetPdgCode()>10000) {AliDebug(1,Form("A photon as track... Track skipped"));continue;}
170 if(!pParticle->GetPDG()) continue;
172 // to be updated for us!!
174 AliDebug(1,Form("Track %i is a %s with charge %i and momentum %f",
175 iTrackN,pParticle->GetPDG()->GetName(),Int_t(pParticle->GetPDG()->Charge()),pParticle->P()));
176 // if(pParticle->GetMother(0)!=-1) continue; //consider only primaries
177 if(pParticle->GetPDG()->Charge()==0||TMath::Abs(Int_t(pParticle->GetPDG()->Charge()))!=3) continue; //to avoid photons from stack...
178 hnvec[0]=pParticle->P();
179 hnvec[1]=pParticle->GetPDG()->Charge();
181 p0.SetMagThetaPhi(pParticle->P(),pParticle->Theta(),pParticle->Phi()); x0.SetXYZ(pParticle->Vx(),pParticle->Vy(),pParticle->Vz());
182 AliRICHHelix helix(x0,p0,TMath::Sign(1,(Int_t)pParticle->GetPDG()->Charge()),b);
183 Int_t iChamber=helix.RichIntersect(AliRICHParam::Instance());
184 if(!iChamber) continue;// no intersection with RICH found
186 hnvec[2]=helix.Ploc().Theta();
187 hnvec[3]=helix.Ploc().Phi();
188 hnvec[4]=helix.PosPc().X();
189 hnvec[5]=helix.PosPc().Y();
191 Double_t dX,dY,dR,dRmip=9999; //min distance between clusters and track position on PC
192 Int_t iMipId=-1; //index of that min distance cluster
193 for(Int_t iClu=0;iClu<pRich->Clus(iChamber)->GetEntries();iClu++){//clusters loop for intersected chamber
194 AliRICHCluster *pClu=(AliRICHCluster*)pRich->Clus(iChamber)->UncheckedAt(iClu);//get pointer to current cluster
195 pClu->DistXY(helix.PosPc(),dX,dY); dR=TMath::Sqrt(dX*dX+dY*dY);//ditance between current cluster and helix intersection with PC
196 if(dR<dRmip){dRmip=dR; hnvec[6]=pClu->X();hnvec[7]=pClu->Y();hnvec[8]=pClu->Q();
197 iMipId=1000000*iChamber+iClu;}//find cluster nearest to the track
199 }//clusters loop for intersected chamber
202 hnvec[9]=recon.ThetaCerenkov(&helix,pRich->Clus(iChamber),iMipId); //search for mean Cerenkov angle for this track
203 hnvec[10]=iMipId;//on return from ThetaCerenkov() contains number of photon candidates accepted
204 hnvec[11]=(Double_t)iMipId;
205 hnvec[12]=(Double_t)iChamber;
206 hnvec[13]=(Double_t)pParticle->GetPdgCode();
207 if(hn) hn->Fill(hnvec);
208 AliDebug(1,Form("FINAL Theta Cerenkov=%f",hnvec[9]));
209 }//stack particles loop
211 pRich->GetLoader()->UnloadRecPoints();
214 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
215 void AliRICHTracker::CalcProb(Double_t thetaCer,Double_t pmod, Double_t *sigmaPID, Double_t *richPID)
217 // Calculates probability to be a electron-muon-pion-kaon-proton
218 // from the given Cerenkov angle and momentum assuming no initial particle composition
219 // (i.e. apriory probability to be the particle of the given sort is the same for all sorts)
220 Double_t height[AliPID::kSPECIES];Double_t totalHeight=0;
221 Double_t thetaTh[AliPID::kSPECIES];
222 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
224 Double_t mass = AliRICHParam::fgMass[iPart];
225 Double_t refIndex=AliRICHParam::Instance()->IdxC6F14(AliRICHParam::EckovMean());
226 Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(refIndex*pmod);
228 if(cosThetaTh>=1) continue;
229 thetaTh[iPart] = TMath::ACos(cosThetaTh);
230 // Double_t sinThetaThNorm = TMath::Sin(thetaTh)/TMath::Sqrt(1-1/(refIndex*refIndex));
231 // Double_t sigmaThetaTh = (0.014*(1/sinThetaThNorm-1) + 0.0043)*1.25;
232 // height[iPart] = TMath::Gaus(thetaCer,thetaTh,sigmaThetaTh);
233 if(sigmaPID[iPart]>0) height[iPart] = TMath::Gaus(thetaCer,thetaTh[iPart],sigmaPID[iPart],kTRUE);
234 else height[iPart] = 0;
235 totalHeight +=height[iPart];
236 AliDebugClass(1,Form(" Particle %s with mass %f with height %f and thetaTH %f",AliPID::ParticleName(iPart),mass,height[iPart],thetaTh[iPart]));
237 AliDebugClass(1,Form(" partial height %15.14f total height %15.14f",height[iPart],totalHeight));
239 if(totalHeight<1e-5) {for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)richPID[iPart]=1.0/AliPID::kSPECIES;return;}
240 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) richPID[iPart] = height[iPart]/totalHeight;
241 Int_t iPartNear = TMath::LocMax(AliPID::kSPECIES,richPID);
242 if(TMath::Abs(thetaCer-thetaTh[iPartNear])>5.*sigmaPID[iPartNear]) for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)richPID[iPart]=1.0/AliPID::kSPECIES;
243 //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
246 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
247 void AliRICHTracker::EsdPrint()
249 // Reads a set of ESD files and print out some information
250 // Arguments: probCut - cut on probability
253 TFile *pFile=TFile::Open("AliESDs.root","read"); if(!pFile) {Printf("ERROR: AliESDs.root does not exist!");return;}
254 TTree *pTr=(TTree*)pFile->Get("esdTree"); if(!pTr) {Printf("ERROR: AliESDs.root, no ESD tree inside!");return;}
255 AliESD *pEsd=new AliESD; pTr->SetBranchAddress("ESD", &pEsd);
257 Int_t iNevt=pTr->GetEntries(); Printf("This ESD contains %i events",iNevt);
258 for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//ESD events loop
260 Int_t iNtracks=pEsd->GetNumberOfTracks(); Printf("ESD contains %i tracks created in Bz=%.2f Tesla",iNtracks,pEsd->GetMagneticField()/10.);
261 for(Int_t iTrk=0;iTrk<iNtracks;iTrk++){//ESD tracks loop
262 AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
263 Double_t dx,dy; pTrack->GetRICHdxdy(dx,dy);
264 Double_t theta,phi; pTrack->GetRICHthetaPhi(theta,phi);
265 Printf("Track %2i Q=%4.1f P=%.3f GeV RICH: ChamberCluster %7i Track-Mip=(%7.2f,%7.2f)=%5.2f cm ThetaCer %7.1f rad",iTrk,pTrack->GetSign(),pTrack->GetP(),
266 pTrack->GetRICHcluster(),dx,dy,TMath::Sqrt(dx*dx+dy*dy),pTrack->GetRICHsignal());
269 delete pEsd; pFile->Close();//close AliESDs.root
271 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
272 void AliRICHTracker::MatrixPrint(Double_t probCut)
274 // 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.
275 // Normally it implies that those 3 ESDs contain only particles of the same sort namly pions, kaons and protons in that order.
276 // Arguments: probCut - cut on probability
278 for(Int_t iFile=0;iFile<3;iFile++){
279 TFile *pFile=TFile::Open(Form("Esd%1i.root",iFile+1),"read"); if(!pFile) {Printf("ERROR: Esd%1i.root does not exist!",iFile+1);return;}
280 TTree *pTr=(TTree*)pFile->Get("esdTree"); if(!pTr) {Printf("ERROR: Esd%1i.root, no ESD tree inside!",iFile+1);return;}
281 AliESD *pEsd=new AliESD; pTr->SetBranchAddress("ESD", &pEsd);
282 Int_t iProtCnt=0,iKaonCnt=0,iPionCnt=0,iUnreconCnt=0,iTrkCnt=0; //counters
284 for(Int_t iEvt=0;iEvt<pTr->GetEntries();iEvt++){//ESD events loop
286 iTrkCnt+=pEsd->GetNumberOfTracks();
287 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//ESD tracks loop
288 AliESDtrack *pTrack = pEsd->GetTrack(iTrk);// get next reconstructed track
289 Double_t dx,dy; pTrack->GetRICHdxdy(dx,dy);
290 Double_t theta,phi; pTrack->GetRICHthetaPhi(theta,phi);
291 Double_t prob[5]; pTrack->GetRICHpid(prob);
292 if(pTrack->GetRICHsignal()>0){
293 if(prob[4]>probCut) iProtCnt++;
294 if(prob[3]>probCut) iKaonCnt++;
295 if((prob[0]+prob[1]+prob[2])>probCut) iPionCnt++;
301 Printf("Bz=%5.2f Events=%i Total tracks=%i No recognized tracks=%i Pion=%i Kaon=%i Proton=%i ProbCut=%.2f",
302 0.1*pEsd->GetMagneticField(),pTr->GetEntries(),iTrkCnt,iUnreconCnt,iPionCnt,iKaonCnt,iProtCnt,probCut);
303 delete pEsd; pFile->Close();//close AliESDs.root
306 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++