]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliDxHFECorrelation.cxx
Merge branch 'feature-movesplit'
[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   Double_t MaxNofEvents = cuts->GetMaxNEventsInPool();
182   Double_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   Double_t MaxNofTracks = (MaxNofEvents+1)*MinNofTracks;
218   Double_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("mc") ||
261         argument.BeginsWith("use-mc")) {
262       fUseMC=true;
263       continue;
264     }
265     if (argument.BeginsWith("system=")) {
266       argument.ReplaceAll("system=", "");
267       if (argument.CompareTo("pp")==0) fSystem=0;
268       else if (argument.CompareTo("Pb-Pb")==0) fSystem=1;
269       else {
270         AliWarning(Form("can not set collision system, unknown parameter '%s'", argument.Data()));
271         // TODO: check what makes sense
272         fSystem=0;
273       }
274       continue;
275     }
276     if (argument.BeginsWith("trigger=")){
277       argument.ReplaceAll("trigger=", "");
278       if (argument.CompareTo("D")==0) { fTriggerParticleType=kD; AliInfo("Trigger on D"); }
279       if (argument.CompareTo("D0")==0) { fTriggerParticleType=kD; AliInfo("Trigger on D");}
280       else if (argument.CompareTo("electron")==0){ fTriggerParticleType=kElectron; AliInfo("trigger on electron");}
281       continue;
282     }   
283     AliWarning(Form("unknown argument '%s'", argument.Data()));
284       
285   }
286
287   return 0;
288 }
289
290 THnSparse* AliDxHFECorrelation::DefineTHnSparse()
291 {
292   //
293   //Defines the THnSparse. For now, only calls CreateControlTHnSparse
294   AliDebug(1, "Creating Corr THnSparse");
295   // here is the only place to change the dimension
296   static const int sizeEventdphi = 7;  
297   InitTHnSparseArray(sizeEventdphi);
298   const double pi=TMath::Pi();
299
300   //TODO: add phi for electron??
301   //                                            0           1       2      3         4     5      6   
302   //                                          D0invmass   PtD0    PhiD0  PtbinD0    Pte   dphi   deta   
303   int         binsEventdphi[sizeEventdphi] = {   200,      1000,   100,    21,     1000,   100,    100};
304   double      minEventdphi [sizeEventdphi] = { 1.5648,      0,       0,     0,       0,  fMinPhi, -2};
305   double      maxEventdphi [sizeEventdphi] = { 2.1648,     100,   2*pi,    20,      100, fMaxPhi, 2};
306   const char* nameEventdphi[sizeEventdphi] = {
307     "D0InvMass",
308     "PtD0",
309     "PhiD0",
310     "PtBinD0",
311     "PtEl",
312     "#Delta#Phi", 
313     "#Delta#eta"
314   };
315
316   TString name;
317   name.Form("%s info", GetName());
318
319
320   return CreateControlTHnSparse(name,sizeEventdphi,binsEventdphi,minEventdphi,maxEventdphi,nameEventdphi);
321
322 }
323
324 THnSparse* AliDxHFECorrelation::CreateControlTHnSparse(const char* name,
325                                                              int thnSize,
326                                                              int* thnBins,
327                                                              double* thnMin,
328                                                              double* thnMax,
329                                                              const char** binLabels) const
330 {
331   //
332   // Creates THnSparse.
333   //
334
335   AliInfo("Setting up THnSparse");
336
337   std::auto_ptr<THnSparseD> th(new THnSparseD(name, name, thnSize, thnBins, thnMin, thnMax));
338   if (th.get()==NULL) {
339     return NULL;
340   }
341   for (int iLabel=0; iLabel<thnSize; iLabel++) {
342     th->GetAxis(iLabel)->SetTitle(binLabels[iLabel]);    
343    
344   }
345   return th.release();
346
347 }
348
349 int AliDxHFECorrelation::AddControlObject(TObject* pObj)
350 {
351   AliInfo("Adding object");
352   /// add control object to list, the base class becomes owner of the object
353   if (!pObj) return -EINVAL;
354   if (!fControlObjects) {
355     fControlObjects=new TList;
356     if (!fControlObjects) return -ENOMEM;
357     fControlObjects->SetOwner();
358   }
359   if (fControlObjects->FindObject(pObj->GetName())) {
360     AliError(Form("ignoring duplicate object '%s' of type %s", pObj->GetName(), pObj->ClassName()));
361     return -EEXIST;
362   }
363   fControlObjects->Add(pObj);
364   return 0;
365 }
366
367 int AliDxHFECorrelation::HistogramEventProperties(int bin)
368 {
369   /// histogram event properties
370   if (!fhEventControlCorr) return 0;
371   fhEventControlCorr->Fill(bin);
372
373   return 0;
374 }
375
376 int AliDxHFECorrelation::Fill(const TObjArray* triggerCandidates, const TObjArray* associatedTracks, const AliVEvent* pEvent)
377 {
378   //
379   // will use AliHFCorrelator to process D0-electron pair and then fill THnSparse.
380   //
381   if (!triggerCandidates || !associatedTracks) return -EINVAL;
382   if (!fControlObjects) {
383     Init();
384   }
385   if (!fControlObjects) {
386     AliError("Initialisation failed, can not fill THnSparse");
387     return -ENODEV;
388   }
389   // set the event to be processed
390   // TODO: change the correlator class to take the const pointer
391   if (!fCorrelator) {
392     AliError("working class instance fCorrelator missing");
393     return -ENODEV;
394   }
395   fCorrelator->SetAODEvent(dynamic_cast<AliAODEvent*>(const_cast<AliVEvent*>(pEvent))); 
396
397   Bool_t correlatorON = fCorrelator->Initialize(); //define the pool for mixing
398   if(!correlatorON) {
399     AliError("AliHFCorrelator didn't initialize the pool correctly or processed a bad event");
400     return 1;
401   }
402
403   TIter itrigger(triggerCandidates);
404   TObject* otrigger=NULL;
405   int ctrigger=-1;
406
407   // For the moment this is very specific to D0-electron correlation. Should be 
408   // changed to be less specific. 
409   while ((otrigger=itrigger())!=NULL) {
410     // loop over trigger D0 particle
411     ctrigger++;
412     AliReducedParticle* ptrigger=dynamic_cast<AliReducedParticle*>(otrigger);
413     if (!ptrigger)  continue;
414
415     Double_t phiTrigger = ptrigger->Phi();
416     Double_t ptTrigger = ptrigger->Pt();
417     Double_t etaTrigger = ptrigger->Eta();
418
419     // set the phi of the D meson in the correct range
420     // TODO: Is this correct to do this??
421     phiTrigger = fCorrelator->SetCorrectPhiRange(phiTrigger);
422     // pass to the object the necessary trigger part parameters
423     fCorrelator->SetTriggerParticleProperties(ptTrigger,phiTrigger,etaTrigger); 
424
425     Bool_t execPool = fCorrelator->ProcessEventPool();
426     if(fUseEventMixing && !execPool) {
427       AliDebug(1,"Mixed event analysis: pool is not ready");
428       continue;
429     }
430     Int_t NofEventsinPool = 1;
431     if(fUseEventMixing) NofEventsinPool = fCorrelator->GetNofEventsInPool();
432                 
433     // loop on events in the pool; if it is SE analysis, stops at one
434     for (Int_t jMix =0; jMix < NofEventsinPool; jMix++){
435       Bool_t analyzetracks = fCorrelator->ProcessAssociatedTracks(jMix, associatedTracks);
436                         
437       if(!analyzetracks) {
438         AliError("AliHFCorrelator::Cannot process the track array");
439         continue;
440       }
441
442       Int_t NofTracks = fCorrelator->GetNofTracks();
443
444       // looping on track candidates
445       for(Int_t iTrack = 0; iTrack<NofTracks; iTrack++){ 
446         Bool_t runcorrelation = fCorrelator->Correlate(iTrack);
447         if(!runcorrelation) continue;
448                         
449         fDeltaPhi = fCorrelator->GetDeltaPhi();
450         fDeltaEta = fCorrelator->GetDeltaEta();
451         
452         AliReducedParticle *assoc = fCorrelator->GetAssociatedParticle();
453
454         FillParticleProperties(ptrigger,assoc,ParticleProperties(),GetDimTHnSparse());
455         fCorrProperties->Fill(ParticleProperties());
456
457       } // loop over electron tracks in event
458     } // loop over events in pool
459   } // loop over trigger particle
460
461   Bool_t updated = fCorrelator->PoolUpdate(associatedTracks);
462   if(fUseEventMixing){
463     if(!updated) AliDebug(1,"Pool was not updated");
464     else {
465       EventMixingChecks(pEvent);
466       AliDebug(1,"Pool was updated");
467     }
468   }
469   return 0;
470 }
471
472
473
474 int AliDxHFECorrelation::FillParticleProperties(AliVParticle* tr, AliVParticle *as, Double_t* data, int dimension) const
475 {
476   // fill the data array from the particle data
477   if (!data) return -EINVAL;
478   AliReducedParticle *ptrigger=(AliReducedParticle*)tr;
479   AliReducedParticle *assoc=(AliReducedParticle*)as;
480   if (!ptrigger || !assoc) return -ENODATA;
481   int i=0;
482   if (dimension!=GetDimTHnSparse()) {
483     // TODO: think about filling only the available data and throwing a warning
484     return -ENOSPC;
485   }
486   if(fTriggerParticleType==kD){
487     data[i++]=ptrigger->GetInvMass();
488     data[i++]=ptrigger->Pt();
489     data[i++]=ptrigger->Phi();
490     data[i++]=ptrigger->GetPtBin(); 
491     data[i++]=assoc->Pt();
492   } 
493   else{
494     data[i++]=assoc->GetInvMass();
495     data[i++]=assoc->Pt();
496     data[i++]=assoc->Phi();
497     data[i++]=assoc->GetPtBin(); 
498     data[i++]=ptrigger->Pt();
499   }
500   data[i++]=GetDeltaPhi();
501   data[i++]=GetDeltaEta();
502
503   return i;
504 }
505
506 void AliDxHFECorrelation::Clear(Option_t * /*option*/)
507 {
508   /// overloaded from TObject: cleanup
509
510   // nothing to be done so far
511   return TObject::Clear();
512 }
513
514 void AliDxHFECorrelation::Print(Option_t */*option*/) const
515 {
516   /// overloaded from TObject: print info
517   cout << "====================================================================" << endl;
518   TObject::Print();
519   if (fHistograms) {
520     fHistograms->Print();
521   }
522 }
523
524 void AliDxHFECorrelation::Draw(Option_t */*option*/)
525 {
526   /// overloaded from TObject: draw histograms
527 }
528
529 TObject* AliDxHFECorrelation::FindObject(const char *name) const
530 {
531   /// overloaded from TObject: find object by name
532   if (fControlObjects) {
533     return fControlObjects->FindObject(name);
534   }
535   return NULL;
536 }
537
538 TObject* AliDxHFECorrelation::FindObject(const TObject *obj) const
539 {
540   /// overloaded from TObject: find object by pointer
541   if (fControlObjects) {
542     return fControlObjects->FindObject(obj);
543   }
544   return NULL;
545 }
546
547 void AliDxHFECorrelation::SaveAs(const char *filename, Option_t */*option*/) const
548 {
549   /// overloaded from TObject: save to file
550   std::auto_ptr<TFile> output(TFile::Open(filename, "RECREATE"));
551   if (!output.get() || output->IsZombie()) {
552     AliError(Form("can not open file %s from writing", filename));
553     return;
554   }
555   output->cd();
556   if (fControlObjects) fControlObjects->Write();
557   output->Close();
558 }
559
560 AliDxHFECorrelation& AliDxHFECorrelation::operator+=(const AliDxHFECorrelation& other)
561 {
562   /// add histograms from another instance
563   // TODO - need to change this to ThnSparse?
564   if (!fHistograms || !other.fHistograms) return *this;
565   
566   for (int i=0; i<kNofHistograms; i++) {
567     if (fHistograms->At(i)==NULL || other.fHistograms->At(i)==NULL) continue;
568     TH1* target=reinterpret_cast<TH1*>(fHistograms->At(i));
569     TH1* source=reinterpret_cast<TH1*>(other.fHistograms->At(i));
570     if (!target || !source) continue;
571     TString name(fHistograms->At(i)->GetName());
572     if (name.CompareTo(target->GetName())!=0) {
573       AliWarning(Form("skipping incompatible objects at position %d: %s vs %s", i, source->GetName(), target->GetName()));
574       continue;
575     }
576     if (source->IsA()!=target->IsA()) {
577       AliWarning(Form("skipping incompatible classes at position %d: %s vs %s", i, source->ClassName(), target->ClassName()));
578       continue;
579     }
580     target->Add(source);
581   }
582   return *this;
583 }
584
585
586 //____________________________  Run checks on event mixing ___________________________________________________
587 void AliDxHFECorrelation::EventMixingChecks(const AliVEvent* pEvent){
588         
589   AliAODEvent *AOD= (AliAODEvent*)(pEvent);
590   AliCentrality *centralityObj = 0;
591   Int_t multiplicity = -1;
592   Double_t MultipOrCent = -1;
593         
594   // get the pool for event mixing
595   if(!fSystem){ // pp
596     multiplicity = AOD->GetNumberOfTracks();
597     MultipOrCent = multiplicity; // convert from Int_t to Double_t
598   }
599   if(fSystem){ // PbPb          
600     centralityObj = ((AliVAODHeader*)AOD->GetHeader())->GetCentralityP();
601     MultipOrCent = centralityObj->GetCentralityPercentileUnchecked("V0M");
602     AliInfo(Form("Centrality is %f", MultipOrCent));
603   }
604         
605   AliAODVertex *vtx = AOD->GetPrimaryVertex();
606   Double_t zvertex = vtx->GetZ(); // zvertex
607
608   AliEventPool *pool = fCorrelator->GetPool();
609
610   ((TH2D*)fControlObjects->FindObject("NofPoolBinCalls"))->Fill(MultipOrCent,zvertex); // number of calls of pool
611   ((TH2D*)fControlObjects->FindObject("EventProps"))->Fill(MultipOrCent,zvertex); // event properties
612   ((TH3D*)fControlObjects->FindObject("EventsPerPoolBin"))->Fill(MultipOrCent,zvertex,pool->NTracksInPool()); // number of events in the pool
613   ((TH3D*)fControlObjects->FindObject("NofTracksPerPoolBin"))->Fill(MultipOrCent,zvertex,pool->GetCurrentNEvents()); // number of calls of pool
614 }
615