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