]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFECorrelation.cxx
adding checks and debugging information
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliDxHFECorrelation.cxx
1 // $Id$\r
2 \r
3 //**************************************************************************\r
4 //* This file is property of and copyright by the ALICE Project            * \r
5 //* ALICE Experiment at CERN, All rights reserved.                         *\r
6 //*                                                                        *\r
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *\r
8 //*                  Sedat Altinpinar <Sedat.Altinpinar@cern.ch>           *\r
9 //*                  Hege Erdal       <hege.erdal@gmail.com>               *\r
10 //*                                                                        *\r
11 //* Permission to use, copy, modify and distribute this software and its   *\r
12 //* documentation strictly for non-commercial purposes is hereby granted   *\r
13 //* without fee, provided that the above copyright notice appears in all   *\r
14 //* copies and that both the copyright notice and this permission notice   *\r
15 //* appear in the supporting documentation. The authors make no claims     *\r
16 //* about the suitability of this software for any purpose. It is          *\r
17 //* provided "as is" without express or implied warranty.                  *\r
18 //**************************************************************************\r
19 \r
20 /// @file   AliDxHFECorrelation.cxx\r
21 /// @author Sedat Altinpinar, Hege Erdal, Matthias Richter\r
22 /// @date   2012-04-25\r
23 /// @brief  Worker class for D0-HF electron correlation\r
24 ///\r
25 \r
26 #include "AliDxHFECorrelation.h"\r
27 #include "AliVParticle.h"\r
28 #include "AliLog.h"\r
29 //#include "AliAnalysisCuts.h"         // required dependency libANALYSISalice.so\r
30 //#include "AliFlowTrackSimple.h"      // required dependency libPWGflowBase.so\r
31 //#include "AliFlowCandidateTrack.h"   // required dependency libPWGflowTasks.so\r
32 //#include "AliCFContainer.h"          // required dependency libCORRFW.so\r
33 #include "TObjArray.h"\r
34 #include "AliHFCorrelator.h"\r
35 #include "AliAODEvent.h"\r
36 #include "AliAODVertex.h"\r
37 #include "TH1D.h"\r
38 #include "TH2D.h"\r
39 #include "TH3D.h"\r
40 #include "THnSparse.h"\r
41 #include "TMath.h"\r
42 #include "TFile.h"\r
43 #include "TCanvas.h"\r
44 #include "TDatabasePDG.h"\r
45 #include "TLorentzVector.h"\r
46 #include "AliReducedParticle.h"\r
47 #include "AliDxHFEParticleSelection.h"\r
48 #include <iostream>\r
49 #include <cerrno>\r
50 #include <memory>\r
51 \r
52 using namespace std;\r
53 \r
54 ClassImp(AliDxHFECorrelation)\r
55 \r
56 AliDxHFECorrelation::AliDxHFECorrelation(const char* name)\r
57   : TNamed(name?name:"AliDxHFECorrelation", "")\r
58   , fHistograms(NULL)  \r
59   , fControlObjects(NULL)\r
60   , fCorrProperties(NULL)\r
61   , fhEventControlCorr(NULL)\r
62   , fCuts(NULL)\r
63   , fUseMC(kFALSE)\r
64   , fCorrelator(NULL)\r
65   , fUseEventMixing(kFALSE)\r
66   , fSystem(0)\r
67   , fMinPhi(-TMath::Pi()/2)\r
68   , fMaxPhi(3*TMath::Pi()/2)\r
69   , fDeltaPhi(0)\r
70   , fDeltaEta(0)\r
71   , fDimThn(0)\r
72   , fCorrArray(NULL)\r
73   , fEventType(0)\r
74 {\r
75   // default constructor\r
76   // \r
77   //\r
78 \r
79 }\r
80 \r
81 const char* AliDxHFECorrelation::fgkEventControlBinNames[]={\r
82   "nEventsAll",\r
83   "nEventsSelected",\r
84   "nEventsD0",\r
85   "nEventsD0e"\r
86 };\r
87 \r
88 AliDxHFECorrelation::~AliDxHFECorrelation()\r
89 {\r
90   // destructor\r
91   //\r
92   //\r
93   if (fHistograms) delete fHistograms;\r
94   fHistograms=NULL;\r
95 \r
96   // NOTE: fControlObjects owns the object, and they are deleted in the\r
97   // destructor of TList\r
98   if (fControlObjects) delete fControlObjects;\r
99   fControlObjects=NULL;\r
100   fCorrProperties=NULL;\r
101   fhEventControlCorr=NULL;\r
102   if(fCorrelator) delete fCorrelator;\r
103   fCorrelator=NULL;\r
104   if(fCorrArray) delete fCorrArray;\r
105   fCorrArray=NULL;\r
106 \r
107   // NOTE: the external object is deleted elsewhere\r
108   fCuts=NULL;\r
109 }\r
110 \r
111 int AliDxHFECorrelation::Init(const char* arguments)\r
112 {\r
113   //\r
114   // Will initialize thnsparse, histogram and AliHFCorrelator\r
115   //\r
116   AliInfo("Initializing correlation objects");\r
117   ParseArguments(arguments);\r
118 \r
119   //----------------------------------------------\r
120   // Setting up THnSparse \r
121   fCorrProperties=DefineTHnSparse();\r
122   AddControlObject(fCorrProperties);\r
123 \r
124   //----------------------------------------------\r
125   // Histogram for storing event information\r
126 \r
127   std::auto_ptr<TH1D> hEventControl(new TH1D("hEventControlCorr", "hEventControlCorr", 10, 0, 10));\r
128   if (!hEventControl.get()) {\r
129     return -ENOMEM;\r
130   }\r
131   int iLabel=0;\r
132   for (iLabel=0; iLabel<kNEventControlLabels; iLabel++)\r
133     hEventControl->GetXaxis()->SetBinLabel(iLabel+1, fgkEventControlBinNames[iLabel]);\r
134 \r
135   fhEventControlCorr=hEventControl.release();\r
136   AddControlObject(fhEventControlCorr);\r
137 \r
138   //----------------------------------------------\r
139   // AliHFCorrelator for Event Mixing and correlation\r
140   // \r
141   // fCuts is the hadron cut object, fSystem to switch between pp or PbPb\r
142   AliHFAssociatedTrackCuts* cuts=dynamic_cast<AliHFAssociatedTrackCuts*>(fCuts);\r
143   if (!cuts) {\r
144     if (fCuts)\r
145       AliError(Form("cuts object of wrong type %s, required AliHFAssociatedTrackCuts", fCuts->ClassName()));\r
146     else\r
147       AliError("mandatory cuts object missing");\r
148     return -EINVAL;\r
149   }\r
150   if (cuts->GetNCentPoolBins()==0 || cuts->GetCentPoolBins()==NULL ||\r
151       cuts->GetNZvtxPoolBins()==0 || cuts->GetZvtxPoolBins()==NULL) {\r
152     // the bin information is used further downstream so it\r
153     // needs to be available in order to continue\r
154     AliError(Form("inavlid object %s: bin configuration is mandatory", cuts->GetName()));\r
155     cuts->Dump();\r
156     return -EINVAL;\r
157   }\r
158   fCorrelator = new AliHFCorrelator("Correlator", cuts, fSystem); \r
159   fCorrelator->SetDeltaPhiInterval(fMinPhi,fMaxPhi); //Correct Phi Interval\r
160   fCorrelator->SetEventMixing(fUseEventMixing);      // mixing Off/On \r
161   fCorrelator->SetAssociatedParticleType(AliHFCorrelator::kElectron);\r
162   // 0: don't calculate d0; 1: return d0; 2: return d0/d0err\r
163   fCorrelator->SetApplyDisplacementCut(kFALSE); \r
164   fCorrelator->SetUseMC(fUseMC);\r
165   fCorrelator->SetUseReco(kTRUE); // Reco/MCTruth\r
166   Bool_t pooldef = fCorrelator->DefineEventPool();\r
167         \r
168   if(!pooldef) AliInfo("Warning:: Event pool not defined properly");\r
169 \r
170 \r
171     // ============================= EVENT MIXING CHECKS ======================================\r
172   // TODO: Not sure if all 4 histos are needed. Keep for now..  \r
173   Int_t MaxNofEvents = cuts->GetMaxNEventsInPool();\r
174   Int_t MinNofTracks = cuts->GetMinNTracksInPool();\r
175   Int_t NofCentBins = cuts->GetNCentPoolBins();\r
176   const Double_t * CentBins = cuts->GetCentPoolBins();\r
177   const Double_t defaultCentBins[] = {0,100};\r
178   if (NofCentBins==0 || CentBins==NULL) {\r
179     NofCentBins=2;\r
180     CentBins=defaultCentBins;\r
181   }\r
182   Int_t NofZVrtxBins = cuts->GetNZvtxPoolBins();\r
183   const Double_t *ZVrtxBins = cuts->GetZvtxPoolBins();\r
184   const Double_t defaultZVrtxBins[] = {-10,10};\r
185   if (NofZVrtxBins==0 || ZVrtxBins==NULL) {\r
186     NofZVrtxBins=2;\r
187     ZVrtxBins=defaultZVrtxBins;\r
188   }\r
189 \r
190   Int_t nofEventPropBins =0;\r
191 \r
192   if(fSystem) nofEventPropBins = 100; // PbPb centrality\r
193   if(!fSystem) nofEventPropBins = NofCentBins; // pp multiplicity\r
194 \r
195   Double_t minvalue = CentBins[0];\r
196   Double_t maxvalue = CentBins[NofCentBins];\r
197   Double_t Zminvalue = ZVrtxBins[0];\r
198   Double_t Zmaxvalue = ZVrtxBins[NofCentBins];\r
199 \r
200   const Double_t Nevents[]={0,2*MaxNofEvents/10,4*MaxNofEvents/10,6*MaxNofEvents/10,8*MaxNofEvents/10,MaxNofEvents};\r
201   const Double_t * events = Nevents;\r
202 \r
203   TH3D * EventsPerPoolBin = new TH3D("EventsPerPoolBin","Number of events in bin pool",NofCentBins,CentBins,NofZVrtxBins,ZVrtxBins,5,events);\r
204   EventsPerPoolBin->GetXaxis()->SetTitle("Centrality/multiplicity ");\r
205   EventsPerPoolBin->GetYaxis()->SetTitle("Z vertex [cm]");\r
206   EventsPerPoolBin->GetZaxis()->SetTitle("Number of events in pool bin");\r
207   if(fUseEventMixing) AddControlObject(EventsPerPoolBin);\r
208 \r
209   Int_t MaxNofTracks = (MaxNofEvents+1)*MinNofTracks;\r
210   Int_t Diff = MaxNofTracks-MinNofTracks;\r
211 \r
212   Double_t Ntracks[]={MinNofTracks,MinNofTracks+Diff/5,MinNofTracks+2*Diff/5,MinNofTracks+3*Diff/5,MinNofTracks+4*Diff/5,MaxNofTracks};\r
213   Double_t  * trackN = Ntracks;\r
214 \r
215   TH3D * NofTracksPerPoolBin = new TH3D("NofTracksPerPoolBin","Number of tracks in bin pool",NofCentBins,CentBins,NofZVrtxBins,ZVrtxBins,5,trackN);\r
216   NofTracksPerPoolBin->GetXaxis()->SetTitle("Centrality/multiplicity ");\r
217   NofTracksPerPoolBin->GetYaxis()->SetTitle("Z vertex [cm]");\r
218   NofTracksPerPoolBin->GetZaxis()->SetTitle("Number of tracks per bin");\r
219 \r
220   if(fUseEventMixing) AddControlObject(NofTracksPerPoolBin);\r
221 \r
222   TH2D * NofPoolBinCalls = new TH2D("NofPoolBinCalls","Number of tracks in bin pool",NofCentBins,CentBins,NofZVrtxBins,ZVrtxBins);\r
223   NofPoolBinCalls->GetXaxis()->SetTitle("Centrality/multiplicity ");\r
224   NofPoolBinCalls->GetYaxis()->SetTitle("Z vertex [cm]");\r
225   if(fUseEventMixing) AddControlObject(NofPoolBinCalls);\r
226 \r
227   TH2D * EventProps = new TH2D("EventProps","Number of tracks in bin pool",nofEventPropBins,minvalue,maxvalue,100,Zminvalue,Zmaxvalue);\r
228   EventProps->GetXaxis()->SetTitle("Centrality/multiplicity ");\r
229   EventProps->GetYaxis()->SetTitle("Z vertex [cm]");\r
230   if(fUseEventMixing) AddControlObject(EventProps);\r
231 \r
232   return 0;\r
233 }\r
234 \r
235 int AliDxHFECorrelation::ParseArguments(const char* arguments)\r
236 {\r
237   // parse arguments and set internal flags\r
238   TString strArguments(arguments);\r
239   auto_ptr<TObjArray> tokens(strArguments.Tokenize(" "));\r
240   if (!tokens.get()) return -ENOMEM;\r
241 \r
242   TIter next(tokens.get());\r
243   TObject* token;\r
244   while ((token=next())) {\r
245     TString argument=token->GetName();\r
246    \r
247     if (argument.BeginsWith("event-mixing")) {\r
248       fUseEventMixing=true;\r
249       continue;\r
250     }\r
251       \r
252     if (argument.BeginsWith("use-mc")) {\r
253       fUseMC=true;\r
254       continue;\r
255     }\r
256     if (argument.BeginsWith("system=")) {\r
257       argument.ReplaceAll("system=", "");\r
258       if (argument.CompareTo("pp")==0) fSystem=0;\r
259       else if (argument.CompareTo("Pb-Pb")==0) fSystem=1;\r
260       else {\r
261         AliWarning(Form("can not set collision system, unknown parameter '%s'", argument.Data()));\r
262         // TODO: check what makes sense\r
263         fSystem=0;\r
264       }\r
265       continue;\r
266     }\r
267     AliWarning(Form("unknown argument '%s'", argument.Data()));\r
268       \r
269   }\r
270 \r
271   return 0;\r
272 }\r
273 \r
274 THnSparse* AliDxHFECorrelation::DefineTHnSparse()\r
275 {\r
276   //\r
277   //Defines the THnSparse. For now, only calls CreateControlTHnSparse\r
278 \r
279   // here is the only place to change the dimension\r
280   static const int sizeEventdphi = 7;  \r
281   InitTHnSparseArray(sizeEventdphi);\r
282   const double pi=TMath::Pi();\r
283 \r
284   //TODO: add phi for electron??\r
285   //                                            0           1       2      3         4     5      6   \r
286   //                                          D0invmass   PtD0    PhiD0  PtbinD0    Pte   dphi   deta   \r
287   int         binsEventdphi[sizeEventdphi] = {   200,      1000,   100,    21,     1000,   100,    100};\r
288   double      minEventdphi [sizeEventdphi] = { 1.5648,      0,       0,     0,       0,  fMinPhi, -2};\r
289   double      maxEventdphi [sizeEventdphi] = { 2.1648,     100,   2*pi,    20,      100, fMaxPhi, 2};\r
290   const char* nameEventdphi[sizeEventdphi] = {\r
291     "D0InvMass",\r
292     "PtD0",\r
293     "PhiD0",\r
294     "PtBinD0",\r
295     "PtEl",\r
296     "#Delta#Phi", \r
297     "#Delta#eta"\r
298   };\r
299 \r
300   TString name;\r
301   name.Form("%s info", GetName());\r
302 \r
303 \r
304   return CreateControlTHnSparse(name,sizeEventdphi,binsEventdphi,minEventdphi,maxEventdphi,nameEventdphi);\r
305 \r
306 }\r
307 \r
308 THnSparse* AliDxHFECorrelation::CreateControlTHnSparse(const char* name,\r
309                                                              int thnSize,\r
310                                                              int* thnBins,\r
311                                                              double* thnMin,\r
312                                                              double* thnMax,\r
313                                                              const char** binLabels) const\r
314 {\r
315   //\r
316   // Creates THnSparse.\r
317   //\r
318 \r
319   AliInfo("Setting up THnSparse");\r
320 \r
321   std::auto_ptr<THnSparseD> th(new THnSparseD(name, name, thnSize, thnBins, thnMin, thnMax));\r
322   if (th.get()==NULL) {\r
323     return NULL;\r
324   }\r
325   for (int iLabel=0; iLabel<thnSize; iLabel++) {\r
326     th->GetAxis(iLabel)->SetTitle(binLabels[iLabel]);    \r
327    \r
328   }\r
329   return th.release();\r
330 \r
331 }\r
332 \r
333 int AliDxHFECorrelation::AddControlObject(TObject* pObj)\r
334 {\r
335   AliInfo("Adding object");\r
336   /// add control object to list, the base class becomes owner of the object\r
337   if (!pObj) return -EINVAL;\r
338   if (!fControlObjects) {\r
339     fControlObjects=new TList;\r
340     if (!fControlObjects) return -ENOMEM;\r
341     fControlObjects->SetOwner();\r
342   }\r
343   if (fControlObjects->FindObject(pObj->GetName())) {\r
344     AliError(Form("ignoring duplicate object '%s' of type %s", pObj->GetName(), pObj->ClassName()));\r
345     return -EEXIST;\r
346   }\r
347   fControlObjects->Add(pObj);\r
348   return 0;\r
349 }\r
350 \r
351 int AliDxHFECorrelation::HistogramEventProperties(int bin)\r
352 {\r
353   /// histogram event properties\r
354   if (!fhEventControlCorr) return 0;\r
355   fhEventControlCorr->Fill(bin);\r
356 \r
357   return 0;\r
358 }\r
359 \r
360 int AliDxHFECorrelation::Fill(const TObjArray* triggerCandidates, const TObjArray* associatedTracks, const AliVEvent* pEvent)\r
361 {\r
362   //\r
363   // will use AliHFCorrelator to process D0-electron pair and then fill THnSparse.\r
364   //\r
365   if (!triggerCandidates || !associatedTracks) return -EINVAL;\r
366   if (!fControlObjects) {\r
367     Init();\r
368   }\r
369   if (!fControlObjects) {\r
370     AliError("Initialisation failed, can not fill THnSparse");\r
371     return -ENODEV;\r
372   }\r
373   // set the event to be processed\r
374   // TODO: change the correlator class to take the const pointer\r
375   if (!fCorrelator) {\r
376     AliError("working class instance fCorrelator missing");\r
377     return -ENODEV;\r
378   }\r
379   fCorrelator->SetAODEvent(dynamic_cast<AliAODEvent*>(const_cast<AliVEvent*>(pEvent))); \r
380 \r
381   Bool_t correlatorON = fCorrelator->Initialize(); //define the pool for mixing\r
382   if(!correlatorON) {\r
383     AliInfo("AliHFCorrelator didn't initialize the pool correctly or processed a bad event");\r
384     return 1;\r
385   }\r
386 \r
387   TIter itrigger(triggerCandidates);\r
388   TObject* otrigger=NULL;\r
389   int ctrigger=-1;\r
390 \r
391   // For the moment this is very specific to D0-electron correlation. Should be \r
392   // changed to be less specific. \r
393   while ((otrigger=itrigger())!=NULL) {\r
394     // loop over trigger D0 particle\r
395     ctrigger++;\r
396     AliReducedParticle* ptrigger=dynamic_cast<AliReducedParticle*>(otrigger);\r
397     if (!ptrigger)  continue;\r
398 \r
399     Double_t phiTrigger = ptrigger->Phi();\r
400     Double_t ptTrigger = ptrigger->Pt();\r
401     Double_t etaTrigger = ptrigger->Eta();\r
402 \r
403     // set the phi of the D meson in the correct range\r
404     // TODO: Is this correct to do this??\r
405     phiTrigger = fCorrelator->SetCorrectPhiRange(phiTrigger);\r
406     // pass to the object the necessary trigger part parameters\r
407     fCorrelator->SetTriggerParticleProperties(ptTrigger,phiTrigger,etaTrigger); \r
408 \r
409     Bool_t execPool = fCorrelator->ProcessEventPool();\r
410     if(fUseEventMixing && !execPool) {\r
411       AliInfo("Mixed event analysis: pool is not ready");\r
412       continue;\r
413     }\r
414     Int_t NofEventsinPool = 1;\r
415     if(fUseEventMixing) NofEventsinPool = fCorrelator->GetNofEventsInPool();\r
416                 \r
417     // loop on events in the pool; if it is SE analysis, stops at one\r
418     for (Int_t jMix =0; jMix < NofEventsinPool; jMix++){\r
419       Bool_t analyzetracks = fCorrelator->ProcessAssociatedTracks(jMix, associatedTracks);\r
420                         \r
421       if(!analyzetracks) {\r
422         AliInfo("AliHFCorrelator::Cannot process the track array");\r
423         continue;\r
424       }\r
425 \r
426       Int_t NofTracks = fCorrelator->GetNofTracks();\r
427 \r
428       // looping on track candidates\r
429       for(Int_t iTrack = 0; iTrack<NofTracks; iTrack++){ \r
430         Bool_t runcorrelation = fCorrelator->Correlate(iTrack);\r
431         if(!runcorrelation) continue;\r
432                         \r
433         fDeltaPhi = fCorrelator->GetDeltaPhi();\r
434         fDeltaEta = fCorrelator->GetDeltaEta();\r
435         \r
436         AliReducedParticle *assoc = fCorrelator->GetAssociatedParticle();\r
437 \r
438         FillParticleProperties(ptrigger,assoc,ParticleProperties(),GetDimTHnSparse());\r
439         fCorrProperties->Fill(ParticleProperties());\r
440 \r
441       } // loop over electron tracks in event\r
442     } // loop over events in pool\r
443   } // loop over trigger particle\r
444 \r
445   Bool_t updated = fCorrelator->PoolUpdate(associatedTracks);\r
446   if(fUseEventMixing){\r
447     if(!updated) AliInfo("Pool was not updated");\r
448     else {\r
449       EventMixingChecks(pEvent);\r
450       AliInfo("Pool was updated");\r
451     }\r
452   }\r
453   return 0;\r
454 }\r
455 \r
456 \r
457 \r
458 int AliDxHFECorrelation::FillParticleProperties(AliVParticle* tr, AliVParticle *as, Double_t* data, int dimension) const\r
459 {\r
460   // fill the data array from the particle data\r
461   if (!data) return -EINVAL;\r
462   AliReducedParticle *ptrigger=(AliReducedParticle*)tr;\r
463   AliReducedParticle *assoc=(AliReducedParticle*)as;\r
464   if (!ptrigger || !assoc) return -ENODATA;\r
465   int i=0;\r
466   if (dimension!=GetDimTHnSparse()) {\r
467     // TODO: think about filling only the available data and throwing a warning\r
468     return -ENOSPC;\r
469   }\r
470   data[i++]=ptrigger->GetInvMass();\r
471   data[i++]=ptrigger->Pt();\r
472   data[i++]=ptrigger->Phi();\r
473   data[i++]=ptrigger->GetPtBin(); \r
474   data[i++]=assoc->Pt();\r
475   data[i++]=GetDeltaPhi();\r
476   data[i++]=GetDeltaEta();\r
477 \r
478   return i;\r
479 }\r
480 \r
481 void AliDxHFECorrelation::Clear(Option_t * /*option*/)\r
482 {\r
483   /// overloaded from TObject: cleanup\r
484 \r
485   // nothing to be done so far\r
486   return TObject::Clear();\r
487 }\r
488 \r
489 void AliDxHFECorrelation::Print(Option_t */*option*/) const\r
490 {\r
491   /// overloaded from TObject: print info\r
492   cout << "====================================================================" << endl;\r
493   TObject::Print();\r
494   if (fHistograms) {\r
495     fHistograms->Print();\r
496   }\r
497 }\r
498 \r
499 void AliDxHFECorrelation::Draw(Option_t */*option*/)\r
500 {\r
501   /// overloaded from TObject: draw histograms\r
502 }\r
503 \r
504 TObject* AliDxHFECorrelation::FindObject(const char *name) const\r
505 {\r
506   /// overloaded from TObject: find object by name\r
507   if (fControlObjects) {\r
508     return fControlObjects->FindObject(name);\r
509   }\r
510   return NULL;\r
511 }\r
512 \r
513 TObject* AliDxHFECorrelation::FindObject(const TObject *obj) const\r
514 {\r
515   /// overloaded from TObject: find object by pointer\r
516   if (fControlObjects) {\r
517     return fControlObjects->FindObject(obj);\r
518   }\r
519   return NULL;\r
520 }\r
521 \r
522 void AliDxHFECorrelation::SaveAs(const char *filename, Option_t */*option*/) const\r
523 {\r
524   /// overloaded from TObject: save to file\r
525   std::auto_ptr<TFile> output(TFile::Open(filename, "RECREATE"));\r
526   if (!output.get() || output->IsZombie()) {\r
527     AliError(Form("can not open file %s from writing", filename));\r
528     return;\r
529   }\r
530   output->cd();\r
531   if (fControlObjects) fControlObjects->Write();\r
532   output->Close();\r
533 }\r
534 \r
535 AliDxHFECorrelation& AliDxHFECorrelation::operator+=(const AliDxHFECorrelation& other)\r
536 {\r
537   /// add histograms from another instance\r
538   // TODO - need to change this to ThnSparse?\r
539   if (!fHistograms || !other.fHistograms) return *this;\r
540   \r
541   for (int i=0; i<kNofHistograms; i++) {\r
542     if (fHistograms->At(i)==NULL || other.fHistograms->At(i)==NULL) continue;\r
543     TH1* target=reinterpret_cast<TH1*>(fHistograms->At(i));\r
544     TH1* source=reinterpret_cast<TH1*>(other.fHistograms->At(i));\r
545     if (!target || !source) continue;\r
546     TString name(fHistograms->At(i)->GetName());\r
547     if (name.CompareTo(target->GetName())!=0) {\r
548       AliWarning(Form("skipping incompatible objects at position %d: %s vs %s", i, source->GetName(), target->GetName()));\r
549       continue;\r
550     }\r
551     if (source->IsA()!=target->IsA()) {\r
552       AliWarning(Form("skipping incompatible classes at position %d: %s vs %s", i, source->ClassName(), target->ClassName()));\r
553       continue;\r
554     }\r
555     target->Add(source);\r
556   }\r
557   return *this;\r
558 }\r
559 \r
560 \r
561 //____________________________  Run checks on event mixing ___________________________________________________\r
562 void AliDxHFECorrelation::EventMixingChecks(const AliVEvent* pEvent){\r
563         \r
564   AliAODEvent *AOD= (AliAODEvent*)(pEvent);\r
565   AliCentrality *centralityObj = 0;\r
566   Int_t multiplicity = -1;\r
567   Double_t MultipOrCent = -1;\r
568         \r
569   // get the pool for event mixing\r
570   if(!fSystem){ // pp\r
571     multiplicity = AOD->GetNTracks();\r
572     MultipOrCent = multiplicity; // convert from Int_t to Double_t\r
573   }\r
574   if(fSystem){ // PbPb          \r
575     centralityObj = AOD->GetHeader()->GetCentralityP();\r
576     MultipOrCent = centralityObj->GetCentralityPercentileUnchecked("V0M");\r
577     AliInfo(Form("Centrality is %f", MultipOrCent));\r
578   }\r
579         \r
580   AliAODVertex *vtx = AOD->GetPrimaryVertex();\r
581   Double_t zvertex = vtx->GetZ(); // zvertex\r
582 \r
583   AliEventPool *pool = fCorrelator->GetPool();\r
584 \r
585   ((TH2D*)fControlObjects->FindObject("NofPoolBinCalls"))->Fill(MultipOrCent,zvertex); // number of calls of pool\r
586   ((TH2D*)fControlObjects->FindObject("EventProps"))->Fill(MultipOrCent,zvertex); // event properties\r
587   ((TH3D*)fControlObjects->FindObject("EventsPerPoolBin"))->Fill(MultipOrCent,zvertex,pool->NTracksInPool()); // number of events in the pool\r
588   ((TH3D*)fControlObjects->FindObject("NofTracksPerPoolBin"))->Fill(MultipOrCent,zvertex,pool->GetCurrentNEvents()); // number of calls of pool\r
589 }\r
590         \r