]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFEParticleSelection.cxx
71ae886e54aae40096bb0d66c48364645793afe3
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliDxHFEParticleSelection.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   AliDxHFEParticleSelection.cxx
21 /// @author Sedat Altinpinar, Hege Erdal, Matthias Richter
22 /// @date   2012-03-19
23 /// @brief  Base class for particle selection
24 ///
25
26 #include "AliDxHFEParticleSelection.h"
27 #include "AliLog.h"
28 #include "AliVEvent.h"
29 #include "AliVParticle.h"
30 #include "TObjArray.h"
31 #include "TList.h"
32 #include "TMath.h"
33 #include "TH1D.h"
34 #include "THnSparse.h"
35 #include "TFile.h"
36 #include <iostream>
37 #include <cerrno>
38 #include <memory>
39
40 using namespace std;
41
42 /// ROOT macro for the implementation of ROOT specific class methods
43 ClassImp(AliDxHFEParticleSelection)
44
45 AliDxHFEParticleSelection::AliDxHFEParticleSelection(const char* name, const char* opt)
46   : TNamed(name?name:"AliDxHFEParticleSelection", name?name:"AliDxHFEParticleSelection")
47   , fOption(opt)
48   , fSelectedTracks(NULL)
49   , fControlObjects(NULL)
50   , fhEventControl(NULL)
51   , fhTrackControl(NULL)
52   , fUseMC(0)
53   , fVerbosity(0)
54   , fDimThn(-1)
55   , fParticleProperties(NULL)
56 {
57   // constructor
58   // 
59   // 
60   // 
61   // 
62 }
63
64 const char* AliDxHFEParticleSelection::fgkEventControlBinNames[]={
65   "nEventsAll",
66   "nEventsSelected",
67   "nEventsWParticle"
68 };
69
70 const char* AliDxHFEParticleSelection::fgkTrackControlBinNames[]={
71   "nTrackAll",
72   "nTrackSelected",
73 };
74
75 AliDxHFEParticleSelection::~AliDxHFEParticleSelection()
76 {
77   // destructor
78   if (fParticleProperties) delete[] fParticleProperties;
79   fParticleProperties=NULL;
80   if (fSelectedTracks) delete fSelectedTracks;
81   fSelectedTracks=NULL;
82   if (fControlObjects) delete fControlObjects;
83   fControlObjects=NULL;
84   fhEventControl=NULL;
85   fhTrackControl=NULL;
86 }
87
88 int AliDxHFEParticleSelection::Init()
89 {
90   //
91   // Init part sel. Calls InitControlObjects()
92   //
93
94   InitControlObjects();
95   return 0;
96 }
97
98
99 int AliDxHFEParticleSelection::InitControlObjects()
100 {
101   //
102   // init control objects
103   // TODO: Change to private now that have Init()?
104   //
105   if (fVerbosity>0) {
106     AliInfo("Setting up control objects");
107   }
108
109   /// init the control objects, can be overloaded by childs which should
110   /// call AliDxHFEParticleSelection::InitControlObjects() explicitly
111   std::auto_ptr<TH1D> hEventControl(new TH1D("hEventControl", "hEventControl", 10, 0, 10));
112   std::auto_ptr<TH1D> hTrackControl(new TH1D("hTrackControl", "hTrackControl", 10, 0, 10));
113
114   fhEventControl=hEventControl.release();
115   for (int iLabel=0; iLabel<kNEventPropertyLabels; iLabel++)
116     fhEventControl->GetXaxis()->SetBinLabel(iLabel+1, fgkEventControlBinNames[iLabel]);
117   AddControlObject(fhEventControl);
118   fhTrackControl=hTrackControl.release();
119   for (int iLabel=0; iLabel<kNTrackPropertyLabels; iLabel++)
120     fhTrackControl->GetXaxis()->SetBinLabel(iLabel+1, fgkTrackControlBinNames[iLabel]);
121   AddControlObject(fhTrackControl);
122
123   return 0;
124 }
125
126 THnSparse* AliDxHFEParticleSelection::CreateControlTHnSparse(const char* name,
127                                                              int thnSize,
128                                                              int* thnBins,
129                                                              double* thnMin,
130                                                              double* thnMax,
131                                                              const char** binLabels) const
132 {
133   //
134   // Creates THnSparse. 
135   //
136
137   AliInfo("Setting up THnSparse");
138
139   std::auto_ptr<THnSparseF> th(new THnSparseF(name, name, thnSize, thnBins, thnMin, thnMax));
140   if (th.get()==NULL) {
141     return NULL;
142   }
143   for (int iLabel=0; iLabel<thnSize; iLabel++) {
144     th->GetAxis(iLabel)->SetTitle(binLabels[iLabel]);    
145    
146   }
147  return th.release();
148
149 }
150
151 THnSparse* AliDxHFEParticleSelection::DefineTHnSparse() 
152 {
153   //
154   // Defines the THnSparse. For now, only calls CreatControlTHnSparse
155
156   //TODO: Should make it more general. Or maybe one can use this here, and skip in PartSelEl?
157
158   // here is the only place to change the dimension
159   const int thnSize = 3;
160   InitTHnSparseArray(thnSize);
161
162   const double Pi=TMath::Pi();
163   TString name;
164   //                                0    1       2
165   //                                Pt   Phi    Eta
166   int         thnBins [thnSize] = { 1000,  200, 500};
167   double      thnMin  [thnSize] = {    0,    0, -1.};
168   double      thnMax  [thnSize] = {  100, 2*Pi,  1.};
169   const char* thnNames[thnSize] = { "Pt","Phi","Eta"};
170
171   name.Form("%s info", GetName());
172
173   return CreateControlTHnSparse(name,thnSize,thnBins,thnMin,thnMax,thnNames);
174 }
175
176 int AliDxHFEParticleSelection::AddControlObject(TObject* pObj)
177 {
178   /// add control object to list, the base class becomes owner of the object
179   if (!pObj) return -EINVAL;
180   if (!fControlObjects) {
181     fControlObjects=new TList;
182     if (!fControlObjects) return -ENOMEM;
183     fControlObjects->SetOwner();
184   }
185   if (fControlObjects->FindObject(pObj->GetName())) {
186     AliError(Form("ignoring duplicate object '%s' of type %s", pObj->GetName(), pObj->ClassName()));
187     return -EEXIST;
188   }
189   if (GetVerbosity()>0) {
190     AliInfo(Form("Adding object '%s' of type %s",pObj->GetName(),pObj->ClassName()));
191   }
192   fControlObjects->Add(pObj);
193   return 0;
194 }
195
196 int AliDxHFEParticleSelection::HistogramEventProperties(int bin)
197 {
198   /// histogram event properties
199   if (!fControlObjects) return 0;
200
201   // TODO: use enums for the bins of the control histogram
202   // for now: 0=all, 1=events with D0s, 2=events with correlated D0s
203   fhEventControl->Fill(bin);
204   return 0;
205 }
206
207 int AliDxHFEParticleSelection::HistogramParticleProperties(AliVParticle* p, int selected)
208 {
209   /// histogram particle properties
210
211   if (!p) return -EINVAL;
212   if (!fControlObjects) return 0;
213
214   // TODO: use enums for the bins of the control histogram
215   fhTrackControl->Fill(kTrackAll);
216   if (selected) fhTrackControl->Fill(kTrackSel);
217   return 0;
218 }
219
220 int AliDxHFEParticleSelection::FillParticleProperties(AliVParticle* p, Double_t* data, int dimension) const
221 {
222   // fill the data array from the particle data
223   if (!data) return -EINVAL;
224   int i=0;
225   if (dimension!=GetDimTHnSparse()) {
226     // TODO: think about filling only the available data and throwing a warning
227     return -ENOSPC;
228   }
229   data[i++]=p->Pt();
230   data[i++]=p->Phi();
231   data[i++]=p->Eta();
232   return i;
233 }
234
235 TObjArray* AliDxHFEParticleSelection::Select(const AliVEvent* pEvent)
236 {
237   /// create selection from 'Tracks' member of the event,
238   /// array contains only pointers but does not own the objects
239   /// object array needs to be deleted by caller
240   if (!pEvent) return NULL;
241   TObjArray* selectedTracks=new TObjArray;
242   if (!selectedTracks) return NULL;
243   int nofTracks=pEvent->GetNumberOfTracks();
244   for (int itrack=0; itrack<nofTracks; itrack++) {
245     AliVParticle* track=pEvent->GetTrack(itrack);
246     int selectionCode=IsSelected(track,pEvent);
247     HistogramParticleProperties(track, selectionCode);
248     if (selectionCode==0) continue;
249     selectedTracks->Add(track);
250   }
251   return selectedTracks;
252 }
253
254 TObjArray* AliDxHFEParticleSelection::Select(TObjArray* pParticles, const AliVEvent* pEvent)
255 {
256   /// create selection from the array of particles,
257   /// array contains only pointers but does not own the objects
258   /// object array needs to be deleted by caller
259   if (!pParticles) return NULL;
260   TObjArray* selectedTracks=new TObjArray;
261   if (!selectedTracks) return NULL;
262   TIter next(pParticles);
263   TObject* pObj=NULL;
264   while ((pObj=next())) {
265     AliVParticle* track=dynamic_cast<AliVParticle*>(pObj);
266     if (!track) continue;
267     int selectionCode=IsSelected(track, pEvent);
268     HistogramParticleProperties(track, selectionCode);
269     if (selectionCode ==0) continue;
270     selectedTracks->Add(track);
271   }
272   return selectedTracks;
273 }
274
275 int AliDxHFEParticleSelection::CheckAndAdd(AliVParticle* /*p*/)
276 {
277   /// check and add track to internal array
278   /// TODO: check if needed
279   return -ENOSYS;
280 }
281
282 int AliDxHFEParticleSelection::IsSelected(AliVParticle* /*p*/, const AliVEvent* /*e*/)
283 {
284   /// check particle if it passes the selection criteria
285   /// childs can overload, by default all tracks are selected
286   return 1;
287 }
288
289 void AliDxHFEParticleSelection::AliDxHFEParticleSelection::Clear(Option_t * /*option*/)
290 {
291   /// inherited from TObject: cleanup
292 }
293
294 void AliDxHFEParticleSelection::Print(Option_t */*option*/) const
295 {
296   /// inherited from TObject: print info
297   cout << "====================================================================" << endl;
298   TNamed::Print();
299   if (fControlObjects) fControlObjects->Print();
300 }
301  
302 void AliDxHFEParticleSelection::SaveAs(const char* filename, Option_t */*option*/) const
303 {
304   /// inherited from TObject: save selection criteria
305   TString fileoption;
306   // TODO: options recreate
307   fileoption="RECREATE";
308   //else fileoption="UPDATE";
309
310   std::auto_ptr<TFile> output(TFile::Open(filename,fileoption));
311   if (!output.get() || output->IsZombie()) {
312     AliError(Form("can not open file %s from writing", filename));
313     return;
314   }
315   output->cd();
316   if (fControlObjects) fControlObjects->Write();
317   output->Close();
318 }
319
320 void AliDxHFEParticleSelection::Draw(Option_t* /*option*/)
321 {
322   /// inherited from TObject: draw content
323
324   // TODO: implement drawing code
325   // - create canvas objects
326   // - plot internal objects
327   // - optionally save canvases to file
328   //
329   // It might be appropriate to have another Draw function taking a
330   // TList as argument and implementing the actual drawing. If this
331   // function is 'static', it can be used stand-alone also from macros
332 }
333
334 TObject* AliDxHFEParticleSelection::FindObject(const char* name) const
335 {
336   /// inherited from TObject: find object by name
337
338   if (fControlObjects) {
339     return fControlObjects->FindObject(name);
340   }
341   return NULL;
342 }
343
344 TObject* AliDxHFEParticleSelection::FindObject(const TObject* obj) const
345 {
346   /// inherited from TObject: find object by pointer
347   if (fControlObjects) {
348     return fControlObjects->FindObject(obj);
349   }
350   return NULL;
351 }