]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHTracker.cxx
New class for PID constants and methods. Changes in all related code (T.Kuhr)
[u/mrichter/AliRoot.git] / RICH / AliRICHTracker.cxx
1 #include "AliRICHTracker.h"
2 #include <AliESD.h>
3 #include <TVector3.h>
4 #include <TTree.h>
5 #include "AliRICH.h"
6 #include "AliRICHHelix.h"
7 #include <AliMagF.h>
8 #include "AliRICHRecon.h"
9 #include <AliStack.h>
10 #include <TParticle.h>
11 #include <TMath.h>
12 ClassImp(AliRICHTracker)
13 //__________________________________________________________________________________________________
14 Int_t AliRICHTracker::PropagateBack(AliESD *pESD)
15 {
16   //invoked by AliRecontruction for RICH
17   //if ESD doesn't contain tracks, try to reconstruct with particle from STACK 
18   //(this case is just to forsee a standalone RICH simulation
19   TNtupleD *hn=0;
20   AliDebug(1,"Start pattern recognition");
21   if(pESD->GetNumberOfTracks()) RecWithESD(pESD);
22   else
23     RecWithStack(hn);
24   AliDebug(1,"Stop pattern recognition");
25
26   return 0; // error code: 0=no error;
27 } //pure virtual from AliTracker
28 //__________________________________________________________________________________________________
29 void AliRICHTracker::RecWithESD(AliESD *pESD)
30 {
31   //recontruction from ESD
32   //
33   Int_t iNtracks=pESD->GetNumberOfTracks();
34   Double_t b=GetFieldMap()->SolenoidField()/10;// magnetic field in Tesla
35   AliDebug(1,Form("Start with %i tracks in %f Tesla field",iNtracks,b));
36   Double_t xb[3],pb[3];//tmp storage for track parameters
37   TVector3 x0(0,0,0); TVector3 p0(0,0,0);//tmp storage for AliRICHHelix
38   
39   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
40   
41   for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//ESD tracks loop
42     AliESDtrack *pTrack = pESD->GetTrack(iTrackN);// get next reconstructed track
43 //  if((pTrack->GetStatus()&AliESDtrack::kTOFout)==0) continue; //ignore tracks not recontructed by TOF
44     pTrack->GetXYZ(xb); 
45     pTrack->GetPxPyPz(pb); 
46     Int_t status=pTrack->GetStatus()&AliESDtrack::kTOFout;//get running track parameters
47     AliDebug(1,Form("Track %i pmod=%f mass=%f stat=%i",iTrackN,pTrack->GetP(),pTrack->GetMass(),status));
48     x0.SetXYZ(xb[0],xb[1],xb[2]); p0.SetXYZ(xb[0],xb[1],xb[2]);
49     AliRICHHelix helix(x0,p0,pTrack->GetSign(),b);   
50     Int_t iChamber=helix.RichIntersect(pRich->P());        
51     AliDebug(1,Form("intersection with %i chamber found",iChamber));
52     if(!iChamber) continue;//intersection with no chamber found
53     
54     Double_t distMip=9999; //min distance between clusters and track position on PC 
55     Int_t iMipId=0; //index of that min distance cluster 
56     for(Int_t iClusN=0;iClusN<pRich->Clusters(iChamber)->GetEntries();iClusN++){//clusters loop for intersected chamber
57       AliRICHcluster *pClus=(AliRICHcluster*)pRich->Clusters(iChamber)->UncheckedAt(iClusN);//get pointer to current cluster
58       Double_t distCurrent=pClus->DistTo(helix.PosPc());//ditance between current cluster and helix intersection with PC
59       if(distCurrent<distMip){distMip=distCurrent;iMipId=iClusN;}//find cluster nearest to the track 
60       
61       AliDebug(1,Form("Ploc (%f,%f,%f) dist= %f",helix.Ploc().Mag(),helix.Ploc().Theta()*TMath::RadToDeg(),
62                                                                     helix.Ploc().Phi()*TMath::RadToDeg(),pClus->DistTo(helix.PosPc())));
63     }////clusters loop for intersected chamber
64     
65     AliDebug(1,Form("Min distance cluster: %i dist is %f",iMipId,distMip));
66     
67     AliRICHRecon recon(&helix,pRich->Clusters(iChamber),iMipId);
68     Double_t thetaCerenkov=recon.ThetaCerenkov(); //search for mean Cerenkov angle for this track
69     AliDebug(1,Form("FINAL Theta Cerenkov=%f",thetaCerenkov));
70     pTrack->SetRICHsignal(thetaCerenkov);
71         
72 //    Double_t richPID[5]={0.2,0.2,0.2,0.2,0.2}; //start with equal probs for (e,mu,pi,k,p)
73 //    CalcProb(thetaCerenkov,p0.Mag(),richPID);
74 //    pTrack->SetRICHpid(richPID);         
75     
76   }//ESD tracks loop
77   AliDebug(1,"Stop.");  
78 } //RecWithESD
79 //__________________________________________________________________________________________________
80 void AliRICHTracker::RecWithStack(TNtupleD *hn)
81 {
82   // reconstruction for particles from STACK
83   //
84   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
85   
86 //  pRich->GetLoader()->GetRunLoader()->LoadHeader();
87   pRich->GetLoader()->GetRunLoader()->LoadKinematics();
88   AliStack *pStack =   pRich->GetLoader()->GetRunLoader()->Stack();
89   Int_t iNtracks=pStack->GetNtrack();
90   AliDebug(1,Form(" Start reconstruction with %i track(s) from Stack",iNtracks));
91   
92   Double_t hnvec[12];
93   
94   Double_t b=GetFieldMap()->SolenoidField()/10;// magnetic field in Tesla
95   AliDebug(1,Form("Start with simulated %i tracks in %f Tesla field",iNtracks,b));
96   TVector3 x0(0,0,0); TVector3 p0(0,0,0);//tmp storage for AliRICHHelix
97   
98
99   if(pRich->GetLoader()->LoadRecPoints()) {AliDebug(1,Form("No clusters found in RICH"));return;}
100   pRich->GetLoader()->TreeR()->GetEntry(0);
101
102   for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//ESD tracks loop
103     TParticle *pParticle = pStack->Particle(iTrackN);
104     AliDebug(1,Form("Track %i is a %s with charge %i and momentum %f",
105             iTrackN,pParticle->GetPDG()->GetName(),Int_t(pParticle->GetPDG()->Charge()),pParticle->P()));
106 //    if(pParticle->GetMother(0)!=-1) continue; //consider only primaries
107     if(pParticle->GetPDG()->Charge()==0||TMath::Abs(Int_t(pParticle->GetPDG()->Charge()))!=3) continue; //to avoid photons from stack...
108     hnvec[0]=pParticle->P();
109     hnvec[1]=pParticle->GetPDG()->Charge();
110     hnvec[2]=pParticle->Theta();
111     hnvec[3]=pParticle->Phi();
112     p0.SetMagThetaPhi(pParticle->P(),pParticle->Theta(),pParticle->Phi());
113     x0.SetXYZ(pParticle->Vx(),pParticle->Vy(),pParticle->Vz());
114     AliRICHHelix helix(x0,p0,TMath::Sign(1,(Int_t)pParticle->GetPDG()->Charge()),b);   
115     Int_t iChamber=helix.RichIntersect(pRich->P());        
116     AliDebug(1,Form("intersection with %i chamber found",iChamber));
117     if(!iChamber) continue;// no intersection with RICH found
118     hnvec[4]=helix.PosPc().X();
119     hnvec[5]=helix.PosPc().Y();
120     Double_t distMip=9999;   //min distance between clusters and track position on PC 
121     Double_t mipX=kBad;      //min distance between clusters and track position on PC 
122     Double_t mipY=kBad;      //min distance between clusters and track position on PC 
123     Double_t chargeMip=kBad; // charge MIP to find
124     Int_t iMipId=kBad;       //index of that min distance cluster 
125     for(Int_t iClusN=0;iClusN<pRich->Clusters(iChamber)->GetEntries();iClusN++){//clusters loop for intersected chamber
126       AliRICHcluster *pClus=(AliRICHcluster*)pRich->Clusters(iChamber)->UncheckedAt(iClusN);//get pointer to current cluster
127       Double_t distCurrent=pClus->DistTo(helix.PosPc());//ditance between current cluster and helix intersection with PC
128       if(distCurrent<distMip){distMip=distCurrent;mipX=pClus->X();
129                                                   mipY=pClus->Y();
130                                                   chargeMip=pClus->Q();iMipId=1000000*iChamber+iClusN;}//find cluster nearest to the track 
131       
132       AliDebug(1,Form("Ploc (%f,%f,%f) dist= %f",helix.Ploc().Mag(),helix.Ploc().Theta()*TMath::RadToDeg(),
133                                                                     helix.Ploc().Phi()*TMath::RadToDeg(),pClus->DistTo(helix.PosPc())));
134     }////clusters loop for intersected chamber
135     
136     AliDebug(1,Form("Min distance cluster: %i dist is %f",iMipId,distMip));
137     hnvec[6]=mipX;hnvec[7]=mipY;
138     hnvec[8]=chargeMip;
139     AliRICHRecon recon(&helix,pRich->Clusters(iChamber),iMipId);
140     Double_t thetaCerenkov=recon.ThetaCerenkov(); //search for mean Cerenkov angle for this track
141     hnvec[9]=thetaCerenkov;
142     hnvec[10]=recon.GetHoughPhotons();
143     hnvec[11]=(Double_t)iMipId;
144     if(hn) hn->Fill(hnvec);
145     AliDebug(1,Form("FINAL Theta Cerenkov=%f",thetaCerenkov));
146 //    pTrack->SetRICHsignal(thetaCerenkov);
147
148 //    Double_t richPID[5]={0.2,0.2,0.2,0.2,0.2}; //start with equal probs for (e,mu,pi,k,p)
149 //    CalcProb(thetaCerenkov,p0.Mag(),richPID);
150     
151   }//ESD tracks loop
152   
153  pRich->GetLoader()->UnloadRecPoints();
154
155   AliDebug(1,"Stop.");  
156 } //RecWithStack
157
158 Int_t AliRICHTracker::LoadClusters(TTree *pTree)
159 {
160 // Load clusters for RICH
161   AliDebug(1,"Start.");  pTree->GetEntry(0);  AliDebug(1,"Stop."); return 0;
162 }
163
164 //__________________________________________________________________________________________________
165 void AliRICHTracker::CalcProb(Double_t thetaCer,Double_t pmod, Double_t *richPID)
166 {
167 // 
168   Double_t height[AliPID::kSPECIES];Double_t totalHeight=0;
169   for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
170     Double_t mass = AliPID::ParticleMass(iPart);
171     Double_t refIndex=AliRICHParam::IndOfRefC6F14(6.755);
172     Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(refIndex*pmod);
173     if(cosThetaTh>=1) {break;}
174     Double_t thetaTh = TMath::ACos(cosThetaTh);
175     Double_t sinThetaThNorm = TMath::Sin(thetaTh)/TMath::Sqrt(1-1/(refIndex*refIndex));
176     Double_t sigmaThetaTh = (0.014*(1/sinThetaThNorm-1) + 0.0043)*1.25;
177     height[iPart] = TMath::Gaus(thetaCer,thetaTh,sigmaThetaTh);
178     totalHeight +=height[iPart];
179   }
180   for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) richPID[iPart] = height[iPart]/totalHeight;    
181 }//CalcProb