]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFEParticleSelectionMCEl.cxx
adding checks and debugging information
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliDxHFEParticleSelectionMCEl.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 //*                  Hege Erdal       <hege.erdal@gmail.com>               *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /// @file   AliDxHFEParticleSelectionMCEl.cxx
20 /// @author Hege Erdal, Matthias Richter
21 /// @date   2012-07-19
22 /// @brief  MC El selection for D0-HFE correlation
23 ///
24
25 #include "AliDxHFEParticleSelectionMCEl.h"
26 #include "AliVParticle.h"
27 #include "AliLog.h"
28 #include "THnSparse.h"
29 #include "AliAODMCParticle.h"
30 #include "TH1F.h"
31 #include "TAxis.h"
32 #include "AliAODTrack.h"
33 #include "AliReducedParticle.h"
34 #include <iostream>
35 #include <cerrno>
36 #include <memory>
37
38 using namespace std;
39 using std::vector;
40
41
42 /// ROOT macro for the implementation of ROOT specific class methods
43 ClassImp(AliDxHFEParticleSelectionMCEl)
44
45 AliDxHFEParticleSelectionMCEl::AliDxHFEParticleSelectionMCEl(const char* opt)
46   : AliDxHFEParticleSelectionEl(opt)
47   , fMCTools()
48   , fPDGnotMCElectron(NULL)
49   , fOriginMother(0)
50   , fResultMC(0)
51 {
52   // constructor
53   // 
54   // 
55   // 
56   // 
57
58   fMCTools.~AliDxHFEToolsMC();
59
60   // TODO: argument scan, build tool options accordingly
61   // e.g. set mc mode first/last, skip control histograms
62   TString toolopt("pdg=11 mc-last");
63   new (&fMCTools) AliDxHFEToolsMC(toolopt);
64 }
65
66 //at the moment not used, keep for now (to be used when plotting)
67 const char* AliDxHFEParticleSelectionMCEl::fgkPDGMotherBinLabels[]={
68   "d",
69   "u",
70   "s",
71   "c",
72   "b",
73   "gluon",
74   "gamma",
75   "#pi^{0}",
76   "#eta",
77   "proton",
78   "others"
79 };
80
81                              
82 const char*  AliDxHFEParticleSelectionMCEl::fgkPDGBinLabels[]={
83   "positron",
84   "electron",
85   "#mu+",
86   "#mu-",
87   "#pi+",
88   "#pi-",
89   "K+",
90   "K-",
91   "proton",
92   "antiproton",
93   "others"
94 };
95
96 AliDxHFEParticleSelectionMCEl::~AliDxHFEParticleSelectionMCEl()
97 {
98   // destructor
99
100   if(fPDGnotMCElectron){
101     delete fPDGnotMCElectron;
102     fPDGnotMCElectron=NULL;
103   }
104
105 }
106
107 int AliDxHFEParticleSelectionMCEl::Init()
108 {
109   //
110   // init function
111   // 
112   int iResult=0;
113   iResult=AliDxHFEParticleSelectionEl::Init();
114   if (iResult<0) return iResult;
115
116   // Histo containing PDG of track which was not MC truth electron
117   fPDGnotMCElectron= new TH1F("fPDGnotMCElectron","PDG of track not MC truth electron",AliDxHFEToolsMC::kNofPDGLabels,-0.5,AliDxHFEToolsMC::kNofPDGLabels-0.5);
118   for (int iLabel=0; iLabel<AliDxHFEToolsMC::kNofPDGLabels; iLabel++)
119     fPDGnotMCElectron->GetXaxis()->SetBinLabel(iLabel+1, fgkPDGBinLabels[iLabel]);
120   AddControlObject(fPDGnotMCElectron);
121   return 0;
122 }
123
124 THnSparse* AliDxHFEParticleSelectionMCEl::DefineTHnSparse()
125 {
126   //
127   // Defines the THnSparse. 
128
129   const int thnSize = 4;
130   InitTHnSparseArray(thnSize);
131   const double Pi=TMath::Pi();
132   TString name;
133   name.Form("%s info", GetName());
134
135   //                           0    1      2     3     
136   //                           Pt   Phi   Eta   mother 
137   int    thnBins[thnSize] = { 1000,  200, 500,    14  };
138   double thnMin [thnSize] = {    0,    0, -1.,  -1.5  };
139   double thnMax [thnSize] = {  100, 2*Pi,  1.,  12.5  };
140   const char* thnNames[thnSize]={
141     "Pt",
142     "Phi",
143     "Eta", 
144     "Mother", //bin==-1: Not MC truth electron
145   };
146  
147   return CreateControlTHnSparse(name,thnSize,thnBins,thnMin,thnMax,thnNames);
148 }
149
150 int AliDxHFEParticleSelectionMCEl::FillParticleProperties(AliVParticle* p, Double_t* data, int dimension) const
151 {
152   // fill the data array from the particle data
153   if (!data) return -EINVAL;
154   AliAODTrack *track=(AliAODTrack*)p;
155   if (!track) return -ENODATA;
156   int i=0;
157   if (dimension!=GetDimTHnSparse()) {
158     // TODO: think about filling only the available data and throwing a warning
159     return -ENOSPC;
160   }
161   data[i++]=track->Pt();
162   data[i++]=track->Phi();
163   data[i++]=track->Eta();
164   data[i++]=fOriginMother; 
165   
166   return i;
167 }
168
169 int AliDxHFEParticleSelectionMCEl::IsSelected(AliVParticle* p, const AliVEvent* pEvent)
170 {
171   /// overloaded from AliDxHFEParticleSelection: check particle
172   /// H: Have changed function. Now doing particle selection first, then run MC over 
173   /// selected tracks. Could configure it to be configurable, but not sure if it
174   /// is needed.  
175   /// Result from normal track selection is returned, result from MC is stored in
176   /// THnSparse. 
177
178   int iResult=0;
179   if (!p || !pEvent){
180     return -EINVAL;
181   }
182   // step 1:
183   // optional MC selection before the particle selection
184   if (fMCTools.MCFirst() && (iResult=CheckMC(p, pEvent))==0) {
185     // histograming?
186     return iResult;
187   }
188
189   // step 2 or 1, depending on sequence:
190   // normal particle selection
191   iResult=AliDxHFEParticleSelectionEl::IsSelected(p, pEvent);
192   if (fMCTools.MCFirst() || iResult==0) return iResult;
193
194   // step 2, only executed if MC check is last
195   // optional MC selection after the particle selection
196   // result stored to be filled into THnSparse
197   // TODO: strictly speaken the particles should be rejected
198   // if not mc selected, however skip this for the moment, because of
199   // the logic outside
200   fResultMC=CheckMC(p, pEvent);
201
202   return iResult;
203 }
204
205 int AliDxHFEParticleSelectionMCEl::CheckMC(AliVParticle* p, const AliVEvent* pEvent)
206 {
207   /// check if MC criteria are fulfilled
208   if (!p || !pEvent){
209     return -EINVAL;
210   }
211   fOriginMother=-1;
212   int iResult=0;
213
214   if (!fMCTools.IsInitialized() && (iResult=fMCTools.InitMCParticles(pEvent))<0) {
215     // TODO: message? but has to be filtered in order to avoid message flood
216     return 0; // no meaningful filtering on mc possible
217   }
218  
219   int pdgParticle=-1;
220   if (fMCTools.RejectByPDG(p,false, &pdgParticle)) {
221     // rejected by pdg
222     // TODO: Move this to fMCTools???? Can this be part of the statistics in the MC class?
223     fPDGnotMCElectron->Fill(fMCTools.MapPDGLabel(pdgParticle));
224     return 0;
225   }
226
227   int pdgMother=0;
228   pdgMother=fMCTools.FindMotherPDG(p,AliDxHFEToolsMC::kGetFirstMother);
229
230   // Particles considered HFE background. can be expanded
231   // Should be created only once 
232   // TODO: that needs to be configured once to avoid performance penalty
233   vector<int> motherPDGs;
234   motherPDGs.push_back(AliDxHFEToolsMC::kPDGpi0); 
235   motherPDGs.push_back(AliDxHFEToolsMC::kPDGeta);
236   motherPDGs.push_back(AliDxHFEToolsMC::kPDGgamma);
237   motherPDGs.push_back(AliDxHFEToolsMC::kPDGJpsi);
238
239   if(fMCTools.RejectByPDG(pdgMother,motherPDGs)){
240     pdgMother=fMCTools.FindMotherPDG(p,AliDxHFEToolsMC::kGetOriginMother);
241     fOriginMother=fMCTools.GetOriginMother();
242   }
243   else{
244     //TODO: Could this be done in a more elegant way?
245     switch(pdgMother){
246     case(AliDxHFEToolsMC::kPDGpi0): fOriginMother=AliDxHFEToolsMC::kNrOrginMother; break;
247     case(AliDxHFEToolsMC::kPDGeta): fOriginMother=AliDxHFEToolsMC::kNrOrginMother+1; break;
248     case(AliDxHFEToolsMC::kPDGgamma): fOriginMother=AliDxHFEToolsMC::kNrOrginMother+2;break;
249     case(AliDxHFEToolsMC::kPDGJpsi): fOriginMother=AliDxHFEToolsMC::kNrOrginMother+3;break;
250     }
251   }
252
253   /*if (fMCTools.RejectByMotherPDG(p)) {
254     // rejected by pdg of original mother
255     // H: want pdg of origin process to be stored in THnSparse
256     // Not sure if this is needed... Need to check, using AliDxHFEToolsMC, who
257     // first mother are, and also what origin is. Use this info here. 
258     return 0;
259     }*/
260
261   return 1;
262 }
263
264 void AliDxHFEParticleSelectionMCEl::Clear(const char* option)
265 {
266   /// clear internal memory
267   fMCTools.Clear(option);
268 }
269
270 AliVParticle *AliDxHFEParticleSelectionMCEl::CreateParticle(AliVParticle* track)
271 {
272
273   AliReducedParticle *part = new AliReducedParticle(track->Eta(), track->Phi(), track->Pt(),track->Charge(),fOriginMother);
274
275   return part;
276
277 }