]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliHFEpidTOF.cxx
Classes for efficiency corrections (Xaver, Roberta)
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliHFEpidTOF.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 #include <TH2F.h>
17 #include <TList.h>
18 #include <TMath.h>
19
20 #include "AliESDtrack.h"
21 #include "AliVParticle.h"
22 #include "AliPID.h"
23
24 #include "AliHFEpidTOF.h"
25 #include "AliHFEpidBase.h"
26
27
28 ClassImp(AliHFEpidTOF)
29   
30 //___________________________________________________________________
31 AliHFEpidTOF::AliHFEpidTOF(const Char_t *name):
32   AliHFEpidBase(name)
33   , fPID(0x0)
34   , fQAList(0x0)
35 {
36   //
37   // Constructor
38   //
39 }
40 //___________________________________________________________________
41 AliHFEpidTOF::AliHFEpidTOF(const AliHFEpidTOF &c):
42   AliHFEpidBase("")
43   , fPID(0x0)
44   , fQAList(0x0)
45 {  
46   // 
47   // Copy operator
48   //
49
50   c.Copy(*this);
51 }
52 //___________________________________________________________________
53 AliHFEpidTOF &AliHFEpidTOF::operator=(const AliHFEpidTOF &ref){
54   //
55   // Assignment operator
56   //
57
58   if(this != &ref){
59     ref.Copy(*this);
60   }
61
62   return *this;
63 }
64 //___________________________________________________________________
65 AliHFEpidTOF::~AliHFEpidTOF(){
66   //
67   // Destructor
68   //
69   if(fPID) delete fPID;
70   if(fQAList){
71     fQAList->Delete();
72     delete fQAList;
73   }
74 }
75 //___________________________________________________________________
76 void AliHFEpidTOF::Copy(TObject &ref) const {
77   //
78   // Performs the copying of the object
79   //
80   AliHFEpidTOF &target = dynamic_cast<AliHFEpidTOF &>(ref);
81
82   target.fPID = fPID;          
83   target.fQAList = fQAList;
84
85   AliHFEpidBase::Copy(ref);
86 }
87 //___________________________________________________________________
88 Bool_t AliHFEpidTOF::InitializePID(){
89   //
90   // InitializePID: TOF experts have to implement code here
91   //
92   return kTRUE;
93 }
94
95 //___________________________________________________________________
96 Int_t AliHFEpidTOF::IsSelected(AliVParticle *_track)
97 {
98
99   //
100   // as of 22/05/2006 :
101   // returns AliPID based on the ESD TOF PID decision
102   // the ESD PID will be checked and if necessary improved 
103   // in the mean time. Best of luck
104   //
105   // returns 10 (== kUnknown)if PID can not be assigned
106   //
107
108   AliESDtrack *track = dynamic_cast<AliESDtrack*>(_track);
109   Long_t status = 0;
110   status = track->GetStatus(); 
111
112   if(!(status & AliESDtrack::kTOFout)) return AliPID::kUnknown;
113   
114   (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(0.);
115
116   Double_t tItrackL = track->GetIntegratedLength();
117   Double_t tTOFsignal = track->GetTOFsignal();
118   
119   if(tItrackL > 0)
120     (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(1.);
121
122   if(tTOFsignal > 0)
123     (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(2.);
124   
125
126   if(tItrackL <=0 || tTOFsignal <=0) return AliPID::kUnknown;
127
128   (dynamic_cast<TH1F *>(fQAList->At(kHistTOFpidFlags)))->Fill(3.);
129   (dynamic_cast<TH1F *>(fQAList->At(kHistTOFsignal)))->Fill(tTOFsignal/1000.);
130   (dynamic_cast<TH1F *>(fQAList->At(kHistTOFlength)))->Fill(tItrackL);
131   // get the TOF pid probabilities
132   Double_t tESDpid[5] = {0., 0., 0., 0., 0.};
133   Float_t tTOFpid_sum = 0.;
134   // find the largest PID probability
135   track->GetTOFpid(tESDpid);
136   Double_t tMAXpid = 0.;
137   Int_t tMAXindex = -1;
138   for(Int_t i=0; i<5; ++i){
139     tTOFpid_sum += tESDpid[i];
140     if(tESDpid[i] > tMAXpid){
141       tMAXpid = tESDpid[i];
142       tMAXindex = i;
143     }
144   }
145   
146   Double_t p = track->GetOuterParam()->P();
147   Double_t beta = (tItrackL/100.)/(TMath::C()*(tTOFsignal/1e12));
148   
149   if(TMath::Abs(tTOFpid_sum - 1) > 0.01) return AliPID::kUnknown;
150   else{
151     // should be the same as AliPID flags
152     
153     (dynamic_cast<TH2F *>(fQAList->At(kHistTOFpid_0+tMAXindex)))->Fill(beta, p);
154     (dynamic_cast<TH2F *>(fQAList->At(kHistTOFpid_beta_v_P)))->Fill(beta, p);
155     return tMAXindex;
156   }
157 }
158
159 //___________________________________________________________________
160 void AliHFEpidTOF::AddQAhistograms(TList *qaList){
161   //
162   // Create QA histograms for TOF PID
163   //
164
165   fQAList = new TList;
166   fQAList->SetName("fTOFqaHistos");
167   fQAList->AddAt(new TH1F("hTOF_flags", "TOF flags;flags (see code for info);counts", 10, -0.25, 4.75), kHistTOFpidFlags);
168   fQAList->AddAt(new TH2F("fTOFbeta_v_P_no","beta -v- P; beta;momentum [GeV/c]", 120, 0, 1.2, 200, 0, 20), kHistTOFpid_beta_v_P);
169   fQAList->AddAt(new TH1F("hTOF_signal", "TOF signal; TOF signal [ns];counts", 1000, 12, 50), kHistTOFsignal);
170   fQAList->AddAt(new TH1F("hTOF_length", "TOF track length; length [cm];counts", 400, 300, 700), kHistTOFlength);
171   fQAList->AddAt(new TH2F("hTOFpid_electron", "TOF reco electron; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid_0);
172   fQAList->AddAt(new TH2F("hTOFpid_muon", "TOF reco muon; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid_1);
173   fQAList->AddAt(new TH2F("hTOFpid_pion", "TOF reco pion; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid_2);
174   fQAList->AddAt(new TH2F("hTOFpid_kaon", "TOF reco kaon; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid_3);
175   fQAList->AddAt(new TH2F("hTOFpid_proton", "TOF reco proton; beta ; momentum [GeV/c]", 120, 0, 1.2, 200, 0, 5), kHistTOFpid_4);
176
177   qaList->AddLast(fQAList);
178 }