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