]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFEParticleSelectionEl.cxx
Update in D-electon correlation code (Hege, Matthias)
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliDxHFEParticleSelectionEl.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE Project            * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  Sedat Altinpinar <Sedat.Altinpinar@cern.ch>           *
9 //*                  Hege Erdal       <hege.erdal@gmail.com>               *
10 //*                                                                        *
11 //* Permission to use, copy, modify and distribute this software and its   *
12 //* documentation strictly for non-commercial purposes is hereby granted   *
13 //* without fee, provided that the above copyright notice appears in all   *
14 //* copies and that both the copyright notice and this permission notice   *
15 //* appear in the supporting documentation. The authors make no claims     *
16 //* about the suitability of this software for any purpose. It is          *
17 //* provided "as is" without express or implied warranty.                  *
18 //**************************************************************************
19
20 /// @file   AliDxHFEParticleSelectionEl.cxx
21 /// @author Sedat Altinpinar, Hege Erdal, Matthias Richter
22 /// @date   2012-03-19
23 /// @brief  D0 selection for D0-HFE correlation
24 ///
25 #include "AliDxHFEParticleSelectionEl.h"
26 #include "AliVParticle.h"
27 #include "AliVEvent.h"
28 #include "AliPID.h"
29 #include "AliPIDResponse.h"
30 #include "AliHFEcontainer.h"
31 #include "AliHFEpid.h"
32 #include "AliHFEpidBase.h"
33 #include "AliHFEtools.h"
34 #include "AliHFEcuts.h"
35 #include "AliAODTrack.h"
36 #include "AliAnalysisDataSlot.h"
37 #include "AliAnalysisDataContainer.h"
38 #include "AliAnalysisManager.h"
39 #include "AliCFManager.h"
40 #include "THnSparse.h"
41 #include "TH1F.h"
42 #include "TAxis.h"
43 #include "TObjArray.h"
44 #include <iostream>
45 #include <cerrno>
46 #include <memory>
47
48 using namespace std;
49
50 /// ROOT macro for the implementation of ROOT specific class methods
51 ClassImp(AliDxHFEParticleSelectionEl)
52
53 AliDxHFEParticleSelectionEl::AliDxHFEParticleSelectionEl(const char* opt)
54   : AliDxHFEParticleSelection("electron", opt)
55   , fPID(NULL)
56   , fElectronProperties(NULL)
57   , fWhichCut(NULL)  
58   , fCuts(NULL)
59   , fCFM(NULL)
60 {
61   // constructor
62   // 
63   // 
64   // 
65   // 
66 }
67
68 AliDxHFEParticleSelectionEl::~AliDxHFEParticleSelectionEl()
69 {
70   // destructor
71   if (fElectronProperties) {
72     delete fElectronProperties;
73     fElectronProperties=NULL;
74   }
75   if(fCFM){
76     delete fCFM;
77     fCFM=NULL;
78   }
79
80   // NOTE: external objects fPID and fCuts are not deleted here
81   fPID=NULL;
82   fCuts=NULL;
83 }
84  
85 int AliDxHFEParticleSelectionEl::InitControlObjects()
86 {
87   /// init control and monitoring objects
88   AliInfo("Electron THnSparse");
89   const int thnSize = 3;
90   const double Pi=TMath::Pi();
91   TString name;// ="e information";
92   //                           0    1       2
93   //                           Pt   Phi    Eta
94   int    thnBins[thnSize] = { 1000,  200, 500};
95   double thnMin [thnSize] = {    0,    0, -1.};
96   double thnMax [thnSize] = {  100, 2*Pi,  1.};
97
98   name.Form("%s info", GetName());
99   std::auto_ptr<THnSparseF> ElectronProperties(new THnSparseF(name, name, thnSize, thnBins, thnMin, thnMax));
100
101   if (ElectronProperties.get()==NULL) {
102     return -ENOMEM;
103   }
104   int axis=0;
105   ElectronProperties->GetAxis(axis++)->SetTitle("Pt");
106   ElectronProperties->GetAxis(axis++)->SetTitle("Phi");
107   ElectronProperties->GetAxis(axis++)->SetTitle("Eta"); 
108
109   fElectronProperties=ElectronProperties.release();
110
111   AddControlObject(fElectronProperties);
112
113   fWhichCut= new TH1F("fWhichCut","effective cut for a rejected particle",6,-0.5,5.5);
114   AddControlObject(fWhichCut);
115
116   //--------Initialize correction Framework and Cuts
117   // Consider moving this, either to separate function or
118   // add a set function for AliCFManager
119   // Do we need this? Can we just call AliHFEcuts::CheckParticleCuts
120   AliInfo("Setting up CFM");
121   fCFM = new AliCFManager;
122   // the setup of cut objects is done in AliHFEcuts::Initialize
123   // the ids used by this class must be the same, the code might be broken if
124   // the sequence in AliHFEcuts::Initialize is changed
125   const Int_t kNcutSteps = AliHFEcuts::kNcutStepsMCTrack + AliHFEcuts::kNcutStepsRecTrack + AliHFEcuts::kNcutStepsDETrack;
126   // reset pointers in the CF manager
127   fCFM->SetNStepParticle(kNcutSteps);
128   for(Int_t istep = 0; istep < kNcutSteps; istep++) {
129     fCFM->SetParticleCutsList(istep, NULL);
130   }
131   if(!fCuts) {
132     AliWarning("Cuts not available. Default cuts will be used");
133     fCuts = new AliHFEcuts;
134     fCuts->CreateStandardCuts();
135   }
136   // TODO: error handling?
137   fCuts->Initialize(fCFM);
138
139   return 0;
140 }
141
142 int AliDxHFEParticleSelectionEl::HistogramParticleProperties(AliVParticle* p, int selectionCode)
143 {
144   /// histogram particle properties
145   if (!p) return -EINVAL;
146   //if (!fControlObjects) return 0;
147   if(selectionCode==0) return  0;
148
149   AliAODTrack *track=(AliAODTrack*)p;
150   Double_t eProperties[]={track->Pt(),track->Phi(),track->Eta()};
151   if(fElectronProperties) fElectronProperties->Fill(eProperties);
152   
153   return 0;
154 }
155
156 int AliDxHFEParticleSelectionEl::IsSelected(AliVParticle* pEl, const AliVEvent*)
157 {
158   /// select El candidates
159   AliAODTrack *track=(AliAODTrack*)pEl;
160  
161
162   //--------track cut selection-----------------------
163   //Using AliHFECuts:
164   // RecKine: ITSTPC cuts  
165   if(!ProcessCutStep(AliHFEcuts::kStepRecKineITSTPC, track)){
166     fWhichCut->Fill(0);
167     return 0;
168   }
169   
170   // RecPrim
171   if(!ProcessCutStep(AliHFEcuts::kStepRecPrim, track)) {
172     fWhichCut->Fill(1);
173     return 0;
174   }
175   
176   // HFEcuts: ITS layers cuts
177   if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsITS, track)) {
178     fWhichCut->Fill(2);
179     return 0;
180   }
181   
182   // HFE cuts: TOF PID and mismatch flag
183   if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTOF, track)) {
184     fWhichCut->Fill(3);
185     return 0;
186   }
187   
188   // HFE cuts: TPC PID cleanup
189   if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTPC, track)){
190     fWhichCut->Fill(4);
191     return 0;
192   } 
193  
194   // HFEcuts: Nb of tracklets TRD0
195   //if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTRD, track)) continue;
196
197
198   //--------PID selection-----------------------
199   AliHFEpidObject hfetrack;
200   hfetrack.SetAnalysisType(AliHFEpidObject::kAODanalysis);
201   hfetrack.SetRecTrack(track);
202
203   // TODO: configurable colliding system
204   //if(IsPbPb()) hfetrack.SetPbPb();
205   hfetrack.SetPP();
206
207   if(fPID && fPID->IsSelected(&hfetrack)) {
208     AliDebug(3,"Inside FilldPhi, electron is selected");
209
210     return 1;
211   }
212   else{
213     fWhichCut->Fill(5);
214     return 0;
215   }
216 }
217
218 void AliDxHFEParticleSelectionEl::SetCuts(TObject* cuts, int level)
219 {
220   /// set cut objects
221   if (level==kCutHFE) {
222     fCuts=dynamic_cast<AliHFEcuts*>(cuts);
223     if (!fCuts && cuts) {
224       AliError(Form("Cut object is not of required type AliHFEcuts but %s", cuts->ClassName()));
225     }
226     return;
227   }
228
229   if (level==kCutPID) {
230     fPID=dynamic_cast<AliHFEpid*>(cuts);
231     if (!fPID && cuts) {
232       AliError(Form("cuts object is not of required type AliHFEpid but %s", cuts->ClassName()));
233     }
234     return;
235   }
236 }
237
238 //________________________________________________________________________
239 Bool_t AliDxHFEParticleSelectionEl::ProcessCutStep(Int_t cutStep, AliVParticle *track)
240 {
241   // Check single track cuts for a given cut step
242   const Int_t kMCOffset = AliHFEcuts::kNcutStepsMCTrack;
243   if(!fCFM->CheckParticleCuts(cutStep + kMCOffset, track)) return kFALSE;
244   //if(!fCuts->CheckParticleCuts(cutStep + kMCOffset, track)) return kFALSE;
245   return kTRUE;
246 }