]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFEParticleSelectionMCD0.cxx
Improve error handling (M. Richter)
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliDxHFEParticleSelectionMCD0.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   AliDxHFEParticleSelectionMCD0.cxx
20 /// @author Hege Erdal, Matthias Richter
21 /// @date   2012-07-19
22 /// @brief  MC D0 selection for D0-HFE correlation
23 ///
24
25 #include "AliDxHFEParticleSelectionMCD0.h"
26 #include "AliAODRecoDecayHF2Prong.h"
27 #include "AliAODTrack.h"
28 #include "AliAODMCParticle.h"
29 #include "AliVParticle.h"
30 #include "AliReducedParticle.h"
31 #include "TH1F.h"
32 #include <iostream>
33 #include <cerrno>
34 #include <memory>
35
36 using namespace std;
37
38 /// ROOT macro for the implementation of ROOT specific class methods
39 ClassImp(AliDxHFEParticleSelectionMCD0)
40
41 AliDxHFEParticleSelectionMCD0::AliDxHFEParticleSelectionMCD0(const char* opt)
42   : AliDxHFEParticleSelectionD0(opt)
43   , fMCTools()
44   , fPDGnotMCD0(NULL)
45   , fResultMC(0)
46   , fOriginMother(0)
47   , fUseKine(kFALSE)
48 {
49   // constructor
50   // 
51   // 
52   // 
53   // 
54
55   // TODO: argument scan, pass only relevant arguments to tools
56   fMCTools.~AliDxHFEToolsMC();
57   TString toolopt("pdg=421 mc-last");
58   new (&fMCTools) AliDxHFEToolsMC(toolopt);
59 }
60
61 AliDxHFEParticleSelectionMCD0::~AliDxHFEParticleSelectionMCD0()
62 {
63   // destructor
64 }
65
66 THnSparse* AliDxHFEParticleSelectionMCD0::DefineTHnSparse()
67 {
68   //
69   // Defines the THnSparse.
70
71   // here is the only place to change the dimension
72   const int thnSize2 = 6;
73   InitTHnSparseArray(thnSize2);
74   const double Pi=TMath::Pi();
75   TString name;
76   name.Form("%s info", GetName());
77
78   //                                 0     1      2       3        4     5
79   //                                 Pt   Phi   Ptbin  D0InvMass  Eta   mother 
80   int         thnBins [thnSize2] = {1000, 200,   15,     200,     500,    10  };
81   double      thnMin  [thnSize2] = {   0,  0,     0,    1.5648,   -1.,  -1.5  };
82   double      thnMax  [thnSize2] = { 100, 2*Pi,  14,    2.1648,    1.,   8.5  };
83   const char* thnNames[thnSize2] = {
84     "Pt",
85     "Phi",
86     "Ptbin", 
87     "D0InvMass", 
88     "Eta",
89     "Mother of D0"  // Bin -1 = not MC truth D0, rest OK
90   };
91
92   // Add Histo displaying pdg of D0 candidates not passing MatchToMC()
93   fPDGnotMCD0= new TH1F("fPDGnotMCD0","PDG of track not MC truth D0",1002,-2.5,999.5);
94   AddControlObject(fPDGnotMCD0);
95
96   return CreateControlTHnSparse(name,thnSize2,thnBins,thnMin,thnMax,thnNames);
97 }
98
99 int AliDxHFEParticleSelectionMCD0::FillParticleProperties(AliVParticle* p, Double_t* data, int dimension) const
100 {
101   // fill the data array from the particle data
102   if (!data) return -EINVAL;
103   AliAODTrack *track=(AliAODTrack*)p;
104   if (!track) return -ENODATA;
105   int i=0;
106   if (dimension!=GetDimTHnSparse()) {
107     // TODO: think about filling only the available data and throwing a warning
108     return -ENOSPC;
109   }
110   data[i++]=track->Pt();
111   data[i++]=track->Phi();
112   data[i++]=AliDxHFEParticleSelectionMCD0::GetPtBin(); 
113   data[i++]=AliDxHFEParticleSelectionMCD0::GetInvMass();
114   data[i++]=track->Eta();
115   data[i++]=fOriginMother; // at the moment not included background. Should expand
116
117   return i;
118 }
119
120 int AliDxHFEParticleSelectionMCD0::IsSelected(AliVParticle* p, const AliVEvent* pEvent)
121 {
122   /// overloaded from AliDxHFEParticleSelection: check particle
123   /// H: Have changed function. Now doing particle selection first, then run MC over 
124   /// selected tracks. Could configure it to be configurable, but not sure if it
125   /// is needed.  
126   /// result from normal track selection is returned, result from MC is stored in
127   /// THnSparse. 
128
129   int iResult=0;
130   fOriginMother=-1;
131
132   // step 1:
133   // MC selection
134   if (fMCTools.MCFirst() && (iResult=CheckMC(p, pEvent))==0) {
135     // histograming?
136     return iResult;
137   }
138
139   // step 2 or 1, depending on sequence:
140   // normal particle selection
141   iResult=AliDxHFEParticleSelectionD0::IsSelected(p, pEvent);
142   if (fMCTools.MCFirst() || iResult==0) return iResult;
143
144   // step 2, only executed if MC check is last
145   // MC selection  - > Should maybe also distinguish between D0 and D0bar
146   // result stored to be filled into THnSparse
147   // TODO: strictly speaken the particles should be rejected
148   // if not mc selected, however skip this for the moment, because of
149   // the logic outside
150   fResultMC=CheckMC(p, pEvent);
151
152   return iResult;
153 }
154
155 int AliDxHFEParticleSelectionMCD0::CheckMC(AliVParticle* p, const AliVEvent* pEvent)
156 {
157   /// check if MC criteria are fulfilled
158   // Check both D0 and D0bar (for now only D0)
159
160   if (!p || !pEvent){
161     return -EINVAL;
162   }
163   int iResult=0;
164
165   if (!fMCTools.IsInitialized() && (iResult=fMCTools.InitMCParticles(pEvent))<0) {
166     // TODO: message? but has to be filtered in order to avoid message flood
167     return 0; // no meaningful filtering on mc possible
168   }
169
170   AliAODRecoDecayHF2Prong *particle = dynamic_cast<AliAODRecoDecayHF2Prong*>(p);
171
172   if(!particle) return 0;
173
174   Int_t pdgDgD0toKpi[2]={AliDxHFEToolsMC::kPDGkaon,AliDxHFEToolsMC::kPDGpion};
175
176   TClonesArray* fMCArray = dynamic_cast<TClonesArray*>(fMCTools.GetMCArray());
177   if(!fMCArray) {cout << "no array" << endl; return -1;}
178
179   // find associated MC particle for D0->Kpi
180   Int_t MClabel=-9999;
181
182   //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx). Checks both D0s and daughters
183   MClabel=particle->MatchToMC(AliDxHFEToolsMC::kPDGD0,fMCArray,2,pdgDgD0toKpi); 
184   
185   if(MClabel<0){
186     // Checking PDG of particle if not MC truth D0
187     // TODO: done the right way??
188     Int_t MCl = p->GetLabel();
189     if(MCl<0) {
190       fPDGnotMCD0->Fill(-2);
191       return 0;
192     }
193     int pdgPart=-1;
194     AliAODMCParticle* aodmcp=0;
195     aodmcp=dynamic_cast<AliAODMCParticle*>(fMCArray->At(MCl));
196     if (aodmcp)
197       pdgPart=TMath::Abs(aodmcp->GetPdgCode());
198     if (pdgPart<0){
199       fPDGnotMCD0->Fill(-1);
200       return 0;
201     }
202     else{
203       fPDGnotMCD0->Fill(pdgPart);
204     }
205     fOriginMother=-1;
206     return 0;
207   }
208
209   fMCTools.SetMClabel(MClabel);
210   fMCTools.FindMotherPDG(p,AliDxHFEToolsMC::kGetOriginMother);
211   fOriginMother=fMCTools.GetOriginMother();
212
213   return 1;
214 }
215
216 void AliDxHFEParticleSelectionMCD0::Clear(const char* option)
217 {
218   /// clear internal memory
219   fMCTools.Clear(option);
220 }
221
222 AliVParticle *AliDxHFEParticleSelectionMCD0::CreateParticle(AliVParticle* track)
223 {
224   //
225   //Created object which contain variables needed for correlation. 
226   //
227
228   AliReducedParticle *part = new AliReducedParticle(track->Eta(), track->Phi(), track->Pt(),AliDxHFEParticleSelectionMCD0::GetInvMass(),AliDxHFEParticleSelectionMCD0::GetPtBin(), fOriginMother);
229
230   return part;
231
232 }