]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFEParticleSelectionD0.cxx
Fix for par file creation on lion with clang (Dario Berzano)
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliDxHFEParticleSelectionD0.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   AliDxHFEParticleSelectionD0.cxx
21 /// @author Sedat Altinpinar, Hege Erdal, Matthias Richter
22 /// @date   2012-03-19
23 /// @brief  D0 selection for D0-HFE correlation
24 ///
25
26 #include "AliDxHFEParticleSelectionD0.h"
27 #include "AliVParticle.h"
28 //#include "AliAnalysisCuts.h"       // required dependency libANALYSISalice.so
29 //#include "AliFlowTrackSimple.h"    // required dependency libPWGflowBase.so
30 //#include "AliFlowCandidateTrack.h" // required dependency libPWGflowTasks.so
31 //#include "AliCFContainer.h"        // required dependency libCORRFW.so
32 #include "AliAODRecoDecayHF2Prong.h" // libPWGHFvertexingHF
33 #include "AliRDHFCutsD0toKpi.h"
34 #include "TObjArray.h"
35 #include "THnSparse.h"
36 #include "TAxis.h"
37 #include "TString.h"
38 #include <iostream>
39 #include <cerrno>
40 #include <memory>
41
42 using namespace std;
43
44 /// ROOT macro for the implementation of ROOT specific class methods
45 ClassImp(AliDxHFEParticleSelectionD0)
46
47 AliDxHFEParticleSelectionD0::AliDxHFEParticleSelectionD0(const char* opt)
48   : AliDxHFEParticleSelection("D0", opt)
49   , fD0Properties(NULL)
50   , fD0Daughter0(NULL)
51   , fD0Daughter1(NULL)
52   , fCuts(NULL)
53   , fFillOnlyD0D0bar(0)
54 {
55   // constructor
56   // 
57   // 
58   // 
59   // 
60   TString strOption(opt);
61   // TODO: one might need a proper argument parsing including
62   // chopping whole string into individual arguments
63   if (strOption.Contains("FillD0D0bar")) fFillOnlyD0D0bar=0;
64   else if (strOption.Contains("FillOnlyD0")) fFillOnlyD0D0bar=1;
65   else if (strOption.Contains("FillOnlyD0bar")) fFillOnlyD0D0bar=2;
66 }
67
68 AliDxHFEParticleSelectionD0::~AliDxHFEParticleSelectionD0()
69 {
70   // destructor
71   if (fD0Properties) {
72     delete fD0Properties;
73     fD0Properties=NULL;
74   }
75   if (fD0Daughter0) {
76     delete fD0Daughter0;
77     fD0Daughter0=NULL;
78   }
79   if (fD0Daughter1) {
80     delete fD0Daughter1;
81     fD0Daughter1=NULL;
82   }
83
84   // Note: external object deleted elsewhere  
85   fCuts=NULL;
86 }
87
88 int AliDxHFEParticleSelectionD0::InitControlObjects()
89 {
90   /// init the control objects, can be overloaded by childs which should
91   /// call AliDxHFEParticleSelection::InitControlObjects() explicitly
92   AliInfo("Setting up THnSparse");
93   TString name;
94   const int thnSize = 4;
95   const double pi=TMath::Pi();
96
97   // TODO: theta?
98   //                             0       1     2      3
99   //                          mass      Pt   Phi    Ptbin
100   int    thnBins[thnSize] = {   200,   1000,  100,   100  };
101   double thnMin [thnSize] = {  1.5648,   0,    0,     0   };
102   double thnMax [thnSize] = {  2.1648, 100,  (2*pi), 100  };
103
104   name.Form("%s info", GetName());
105   std::auto_ptr<THnSparseF> D0Properties(new THnSparseF(name, name, thnSize, thnBins, thnMin, thnMax));
106
107   if (D0Properties.get()==NULL) {
108     return -ENOMEM;
109   }
110   int axis=0;
111   D0Properties->GetAxis(axis++)->SetTitle("D0 inv mass");
112   D0Properties->GetAxis(axis++)->SetTitle("Pt");
113   D0Properties->GetAxis(axis++)->SetTitle("Phi"); 
114   D0Properties->GetAxis(axis++)->SetTitle("Ptbin"); 
115
116   fD0Properties=D0Properties.release();
117   AddControlObject(fD0Properties);
118
119   //Adding control objects for the daughters
120   InitControlObjectsDaughters("pi information",0);
121   InitControlObjectsDaughters("K information",1);
122
123   return AliDxHFEParticleSelection::InitControlObjects();
124 }
125
126 int AliDxHFEParticleSelectionD0::InitControlObjectsDaughters(TString name, int daughter)
127 {
128   //Setting up Control objects for the daughters.
129   AliInfo("Setting up daughter THnSparse");
130
131   const int thnSize2 = 5;
132   const double Pi=TMath::Pi();
133   //                           0    1      2      3          4
134   //                           Pt   Phi   Ptbin  D0InvMass  Eta
135   int    thnBins[thnSize2] = { 1000,  200, 21,     200,     500};
136   double thnMin [thnSize2] = {    0,    0,  0,    1.5648,   -1.};
137   double thnMax [thnSize2] = {  100, 2*Pi, 20,    2.1648,    1.};
138
139   std::auto_ptr<THnSparseF> DaughterProperties(new THnSparseF(name, name, thnSize2, thnBins, thnMin, thnMax));
140
141   if (DaughterProperties.get()==NULL) {
142     return -ENOMEM;
143   }
144   int axis=0;
145   DaughterProperties->GetAxis(axis++)->SetTitle("Pt");
146   DaughterProperties->GetAxis(axis++)->SetTitle("Phi");
147   DaughterProperties->GetAxis(axis++)->SetTitle("Ptbin"); 
148   DaughterProperties->GetAxis(axis++)->SetTitle("D0InvMass"); 
149   DaughterProperties->GetAxis(axis++)->SetTitle("Eta"); 
150
151   if(daughter==0){ 
152     fD0Daughter0=DaughterProperties.release();
153     AddControlObject(fD0Daughter0);
154   }
155   
156   if(daughter==1){
157     fD0Daughter1=DaughterProperties.release();
158     AddControlObject(fD0Daughter1);
159   }
160   return 0;
161 }
162
163 int AliDxHFEParticleSelectionD0::HistogramParticleProperties(AliVParticle* p, int selectionCode)
164 {
165
166   /// histogram particle properties
167   if (!p) return -EINVAL;
168
169   // fill the common histograms
170   AliDxHFEParticleSelection::HistogramParticleProperties(p, selectionCode);
171
172   // no daughters to fill if 0 (= no candidate)
173   if (selectionCode==0) return 0;
174
175   AliAODRecoDecayHF2Prong* part=dynamic_cast<AliAODRecoDecayHF2Prong*>(p);
176
177   if(!part) return 0;
178   // Convention: 1. daughter is postive track, 2. = negative
179   AliAODTrack *prongpos=(AliAODTrack*)part->GetDaughter(0);
180   AliAODTrack *prongneg=(AliAODTrack*)part->GetDaughter(1);
181
182   if(!prongpos || !prongneg) {
183     return 0;
184   }
185  
186   // Only D0s are filled 
187   // TODO: Also include D0bar
188   if ((selectionCode==1 || selectionCode==3) && fFillOnlyD0D0bar<2) {
189     Double_t invmassD0 = part->InvMassD0();
190     AliDebug(3,Form("pt %f,  nr pt bins %d",part->Pt(),fCuts->GetNPtBins()));
191     Int_t ptbin=fCuts->PtBin(part->Pt());
192     Double_t D0Properties[] = {invmassD0,part->Pt(),part->Phi(),ptbin};
193     Double_t KProperties[]={prongneg->Pt(),prongneg->Phi(),ptbin, invmassD0,prongneg->Eta()};
194     Double_t piProperties[]={prongpos->Pt(),prongpos->Phi(),ptbin,invmassD0,prongpos->Eta()};
195
196     if(fD0Properties) fD0Properties->Fill(D0Properties);
197     if(fD0Daughter0) fD0Daughter0->Fill(piProperties);
198     if(fD0Daughter1) fD0Daughter1->Fill(KProperties);
199   }
200
201   return 0;
202 }
203
204 TObjArray* AliDxHFEParticleSelectionD0::Select(TObjArray* pTracks, const AliVEvent *pEvent)
205 {
206   /// create selection, array contains only pointers but does not own the objects
207   /// object array needs to be deleted by caller
208   if (!pTracks) return NULL;
209   TObjArray* selectedTracks=new TObjArray;
210   if (!selectedTracks) return NULL;
211   TIter itrack(pTracks);
212   TObject* pObj=NULL;
213   while ((pObj=itrack())!=NULL) {
214     AliVParticle* track=dynamic_cast<AliVParticle*>(pObj);
215     if (!track) continue;
216     int selectionCode=IsSelected(track,pEvent);
217     HistogramParticleProperties(track, selectionCode);
218     //TODO: Also add selection for D0bar
219
220     // Add track if it is either defined as D0(selectionCode==1) or both 
221     // D0bar and a D0 (selectionCode==3)
222     if (! ((selectionCode==1 || selectionCode==3) && fFillOnlyD0D0bar<2)) continue;
223     selectedTracks->Add(track);
224   }
225   return selectedTracks;
226 }
227
228 int AliDxHFEParticleSelectionD0::IsSelected(AliVParticle* p, const AliVEvent* pEvent)
229 {
230   /// TODO: implement specific selection of D0 candidates
231   /// Could also return values based on where where selection "failed"
232   int selectionCode=0;
233
234   AliAODRecoDecayHF2Prong *d0 = dynamic_cast<AliAODRecoDecayHF2Prong*>(p);
235   if(d0->GetSelectionMap()) if(!d0->HasSelectionBit(AliRDHFCuts::kD0toKpiCuts)){
236       AliDebug(1,"Skip D0 from Dstar");
237       return 0; //skip the D0 from Dstar
238     }
239
240   // TODO: the cuts instance should be const but the function definition of
241   // AliRDHFCuts::IsSelected does not allow this
242   AliRDHFCuts* cuts=const_cast<AliRDHFCuts*>(fCuts);
243   if (!cuts) {
244     selectionCode=1;
245   } else if(cuts->IsInFiducialAcceptance(d0->Pt(),d0->Y(421)) ) {
246     //    if(cuts->IsSelected(d0,AliRDHFCuts::kTracks,pEvent))fNentries->Fill(6);       
247     Int_t ptbin=cuts->PtBin(d0->Pt());
248     if(ptbin==-1) {
249       AliDebug(1,"Pt out of bounds");
250       return 0;
251     } //out of bounds
252
253     // TODO: the aod pointer should also be const but the function definition of
254     // AliRDHFCuts::IsSelected does not allow this
255     AliAODEvent* aod=NULL;
256     if (pEvent) aod=dynamic_cast<AliAODEvent*>(const_cast<AliVEvent*>(pEvent));
257   
258     // Selected. Return 0 (none), 1(D0), 2(D0bar) or 3 (both)
259     selectionCode=cuts->IsSelected(d0,AliRDHFCuts::kAll,aod); 
260
261     AliDebug(1,Form("Candidate is %d \n", selectionCode));
262     TObjArray daughters;
263     daughters.AddAt((AliAODTrack*)d0->GetDaughter(0),0);
264     daughters.AddAt((AliAODTrack*)d0->GetDaughter(1),1);
265
266     //check daughters
267     if(!daughters.UncheckedAt(0) || !daughters.UncheckedAt(1)) {
268       AliDebug(1,"at least one daughter not found!");
269       daughters.Clear();
270       return 0;
271     }
272   }
273
274   return selectionCode;
275 }
276
277 void AliDxHFEParticleSelectionD0::SetCuts(TObject* cuts, int /*level*/)
278 {
279   /// set cuts objects
280   fCuts=dynamic_cast<AliRDHFCuts*>(cuts);
281   if (!fCuts && cuts) {
282     AliError(Form("cuts object is not of required type AliRDHFCuts but %s", cuts->ClassName()));
283   }
284 }