]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/Centrality/AliAnalysisTaskHIMultCorr.cxx
Updates for the centralized automatic QA (Ionut)
[u/mrichter/AliRoot.git] / PWGPP / Centrality / AliAnalysisTaskHIMultCorr.cxx
1 //-*- Mode: C++ -*-
2
3 #if __GNUS__ >= 3
4 using namespace std;
5 #endif
6
7 #include "TFile.h"
8 #include "TChain.h"
9 #include "TTree.h"
10 #include "TH1F.h"
11 #include "TH2F.h"
12 #include "TCanvas.h"
13 #include "TStopwatch.h"
14 #include "TMath.h"
15
16 #include "AliAnalysisTask.h"
17 #include "AliAnalysisManager.h"
18 #include "AliTracker.h" 
19 #include "AliESDEvent.h"
20 #include "AliESDInputHandler.h"
21 #include "AliStack.h"
22 #include "AliMCEvent.h"
23 #include "AliMCEventHandler.h"
24 #include "AliESDtrackCuts.h"
25 #include "AliKineTrackCuts.h"
26 #include "AliMCParticle.h"
27 #include "AliESDVZERO.h"
28 #include "AliAnalysisTaskHIMultCorr.h"
29 #include "AliCentrality.h"
30
31 #include "AliESDVZERO.h"
32 #include "AliTriggerAnalysis.h"
33 #include "AliMultiplicity.h"
34
35 #include "AliMultiplicityCorrelations.h"
36
37 // Task for HI Multiplicity correlation checks
38 // Author: Jochen Thaeder <jochen@thaeder.de>
39
40 ClassImp(AliAnalysisTaskHIMultCorr)
41
42 #define USE_STREAMER 0
43
44 /*
45  * ---------------------------------------------------------------------------------
46  *                            Constructor / Destructor
47  * ---------------------------------------------------------------------------------
48  */
49
50 //________________________________________________________________________
51 AliAnalysisTaskHIMultCorr::AliAnalysisTaskHIMultCorr(const char *name) :
52   AliAnalysisTaskSE(name),
53   fpcstream(NULL),
54   fIsMC(kFALSE),
55   fHStat(NULL),
56   fOutList(NULL),
57   fESD(NULL), fESDTrackCuts(NULL), fESDTrackCuts2(NULL),
58   fUseCentralitySel(0),
59   fCentralityBin(-1),
60   fCentralitySPDBin(-1),
61   fCentralityVZEROBin(-1),
62   fCentralitySPD(-1.),
63   fCentralityVZERO(-1.),
64   fMaxVertexZ(30.),
65   fTriggerAnalysis(NULL),
66   fCorrObj(NULL), 
67   fCorrObjCent0(NULL), fCorrObjCent1(NULL), fCorrObjCent2(NULL) {
68   // Constructor   
69
70 #if USE_STREAMER
71   fpcstream = new TTreeSRedirector("eventInfoCorr.root");
72 #endif
73
74   // Output slot #1 writes into a TList container
75   DefineOutput(1, TList::Class());
76 }
77
78 //________________________________________________________________________
79 AliAnalysisTaskHIMultCorr::~AliAnalysisTaskHIMultCorr() {
80   // Destructor
81
82   if (fTriggerAnalysis)
83     delete fTriggerAnalysis;
84   fTriggerAnalysis = NULL;
85
86 #if USE_STREAMER
87   if (fpcstream)
88     delete fpcstream;
89   fpcstream = NULL;
90 #endif
91
92 }
93
94 /*
95  * ---------------------------------------------------------------------------------
96  *                                    Methods
97  * ---------------------------------------------------------------------------------
98  */
99
100 //________________________________________________________________________
101 void AliAnalysisTaskHIMultCorr::UserCreateOutputObjects() {
102   // Create histograms
103   // Called once
104
105   Bool_t oldStatus = TH1::AddDirectoryStatus();
106   TH1::AddDirectory(kFALSE);
107
108   fOutList = new TList;
109   fOutList->SetName(GetName()) ;
110   fOutList->SetOwner(kTRUE) ;
111
112   fTriggerAnalysis = new AliTriggerAnalysis;
113
114   // ------------------------------------------------------------------
115   // -- Correlations
116   // ------------------------------------------------------------------
117
118   fCorrObj = new AliMultiplicityCorrelations;
119   fCorrObj->Initialize("All");
120   fOutList->Add(fCorrObj->GetHistList());
121
122   fCorrObjCent0 = new AliMultiplicityCorrelations;
123   fCorrObjCent0->Initialize("Cent_0-5");
124   //  fCorrObjCent0->SetCleanSample(1700.,2900.);      // Cleaning of Sample
125   fOutList->Add(fCorrObjCent0->GetHistList());
126
127   fCorrObjCent1 = new AliMultiplicityCorrelations;
128   fCorrObjCent1->Initialize("Cent_80-90");
129   fOutList->Add(fCorrObjCent1->GetHistList());
130   
131   fCorrObjCent2 = new AliMultiplicityCorrelations;
132   fCorrObjCent2->Initialize("Cent_70-80");
133   //  fCorrObjCent2->SetCleanSample(0.,180.);         // Cleaning of Sample
134   fOutList->Add(fCorrObjCent2->GetHistList());
135
136   if (fIsMC) {
137     fCorrObj->SetIsMC();
138     fCorrObjCent0->SetIsMC();
139     fCorrObjCent1->SetIsMC();
140     fCorrObjCent2->SetIsMC();
141   }
142
143   if (fESDTrackCuts) {
144     fCorrObj->SetESDTrackCuts(fESDTrackCuts);
145     fCorrObjCent0->SetESDTrackCuts(fESDTrackCuts);
146     fCorrObjCent1->SetESDTrackCuts(fESDTrackCuts);
147     fCorrObjCent2->SetESDTrackCuts(fESDTrackCuts);
148   }
149   else
150     AliError("No ESD trackCuts!!");
151
152   if (fESDTrackCuts2) {
153     fCorrObj->SetESDTrackCuts2(fESDTrackCuts2);
154     fCorrObjCent0->SetESDTrackCuts2(fESDTrackCuts2);
155     fCorrObjCent1->SetESDTrackCuts2(fESDTrackCuts2);
156     fCorrObjCent2->SetESDTrackCuts2(fESDTrackCuts2);
157   }
158
159   // ------------------------------------------------------------------
160   // -- Centrality Distributions
161   // ------------------------------------------------------------------
162   
163   fOutList->Add(new TH1F("fVZEROCentrality",   "VZERO Centrality Distribution;Centrality;N_{Events}",     102,0.,101.));
164   fOutList->Add(new TH1F("fVZEROCentralityBin","VZERO Centrality Bin Distribution;Centrality;N_{Events}", 102,0.,101.));
165
166   fOutList->Add(new TH1F("fSPDCentrality",   "SPD Centrality Distribution;Centrality;N_{Events}",         102,0.,101.));
167   fOutList->Add(new TH1F("fSPDCentralityBin","SPD Centrality Bin Distribution;Centrality;N_{Events}",     102,0.,101.));
168
169   fOutList->Add(new TH2F("fCorrCentrality",
170                          "Centrality - SPD vs VZERO;Centrality_{SPD};Centrality_{VZERO}",      102,0.,101., 102,0.,101.));
171   fOutList->Add(new TH2F("fCorrCentralityBin",
172                          "Centrality Bin -  SPD vs VZERO;Centrality_{SPD};Centrality_{VZERO}", 102,0.,101., 102,0.,101.));
173
174   // ------------------------------------------------------------------
175   // -- Vertex : TPC - Tracks
176   // ------------------------------------------------------------------
177   
178   fOutList->Add(new TH2F("fCorrVertexTracks",
179                          "Tracks Vertex - Status vs Contributers;Status;N_{contributers}",     2,-1.,3., 500,0.,3000.));
180   fOutList->Add(new TH2F("fCrossCorrVertexTracks",
181                          "Tracks Vertex - TPC Status vs Contributers;Status;N_{contributers}", 2,-1.,3., 500,0.,3000.));
182
183   fOutList->Add(new TH2F("fCorrVertexTPC",
184                          "TPC Vertex - Status vs Contributers;Status;N_{contributers}",        2,-1.,3., 500,0.,3000.));
185   fOutList->Add(new TH2F("fCrossCorrVertexTPC",
186                          "TPC Vertex - Tracks Status vs Contributers;Status;N_{contributers}", 2,-1.,3., 500,0.,3000.));
187
188   fOutList->Add(new TH2F("fCorrVertexNCont",
189                          "N Contributers - Vertex Tracks vs Vertex TPC;N_{contributers}^{Tracks};N_{contributers}^{TPC}", 
190                          500,0.,3000., 500,0.,3000.));
191
192   fOutList->Add(new TH1F("fDeltaVx","Vertex Tracks - Vertex TPC #Deltax;#Delta X [cm]", 201,-10.,10.));
193   fOutList->Add(new TH1F("fDeltaVy","Vertex Tracks - Vertex TPC #Deltay;#Delta Y [cm]", 201,-10.,10.));
194   fOutList->Add(new TH1F("fDeltaVz","Vertex Tracks - Vertex TPC #Deltaz;#Delta z [cm]", 201,-10.,10.));
195
196   // ------------------------------------------------------------------
197   // -- Vertex Z : TPC - Tracks - SPD
198   // ------------------------------------------------------------------
199
200   fOutList->Add(new TH2F("fVzTPCvsSPD",   
201                          "Z - Vertex TPC vs Vertex SPD;Vz_{TPC} [cm];Vz_{SPD} [cm]",       401,-20.,20., 401,-20.,20.));
202   fOutList->Add(new TH2F("fVzTPCvsTracks", 
203                          "Z - Vertex TPC vs Vertex Tracks;Vz_{TPC} [cm];Vz_{Tracks} [cm]", 401,-20.,20., 401,-20.,20.));
204   fOutList->Add(new TH2F("fVzSPDvsTracks",
205                          "Z - Vertex SPD vs Vertex Tracks;Vz_{SPD} [cm];Vz_{Tracks} [cm]", 401,-20.,20., 401,-20.,20.));
206
207   // ------------------------------------------------------------------
208   // -- Cuts 
209   // ------------------------------------------------------------------
210
211   fOutList->Add(new TH1F("fHStat","Cut statistics", 20,0.,19));
212   fHStat = dynamic_cast<TH1F*>(fOutList->Last());
213
214   TH1::AddDirectory(oldStatus);
215 }
216
217 //________________________________________________________________________
218 void AliAnalysisTaskHIMultCorr::UserExec(Option_t *) {
219   // Called for each event
220   
221 #if USE_STREAMER
222     Int_t   i1 = -1;
223     Int_t   i2 = -1;
224     Float_t f1 = -1.;
225     Float_t f2 = -1.;
226     Float_t f3 = -1.;
227 #endif
228
229   // -- Setup Event
230   // ----------------
231   if ( !SetupEvent() ) {
232 #if USE_STREAMER
233     (*fpcstream) <<"event"<<
234       "nTracksTPC=" << i1 <<
235       "nTracks="    << i2 <<
236       "v0A="        << f1 <<
237       "v0C="        << f2 <<
238       "v0Corr="     << f3 << "\n";
239 #endif
240     return;
241   }
242     
243   // -- Centrality dependence
244   // --------------------------
245   if (fCentralityBin < 0) { 
246 #if USE_STREAMER
247     (*fpcstream) <<"event"<<
248       "nTracksTPC=" << i1 <<
249       "nTracks="    << i2 <<
250       "v0A="        << f1 <<
251       "v0C="        << f2 <<
252       "v0Corr="     << f3 << "\n";
253 #endif
254     return;
255   }
256     
257   // -- Process Event
258   // ------------------
259   Int_t iResult = 0;
260   if (fCorrObj) iResult = fCorrObj->ProcessEvent(fESD);
261   if (iResult == -3) {
262 #if USE_STREAMER
263     (*fpcstream) <<"event"<<
264       "nTracksTPC=" << i1 <<
265       "nTracks="    << i2 <<
266       "v0A="        << f1 <<
267       "v0C="        << f2 <<
268       "v0Corr="     << f3 << "\n";
269 #endif
270     return;
271   }
272   
273   // -- After "Jochen's cut"
274   fHStat->Fill(7);
275   fHStat->Fill(17);
276
277   if      (fCentralityBin == 0  && fCorrObjCent0) fCorrObjCent0->ProcessEvent(fESD);
278   else if (fCentralityBin == 80 && fCorrObjCent1) fCorrObjCent1->ProcessEvent(fESD);
279   else if (fCentralityBin == 70 && fCorrObjCent2) fCorrObjCent2->ProcessEvent(fESD);
280
281   // -- Marian's Streamer
282   // ----------------------
283 #if USE_STREAMER
284   if (fCorrObj) {
285     i1 = fCorrObj->GetNTracksTPC();
286     i2 = fCorrObj->GetNTracks();
287     f1 = fCorrObj->GetVZEROA();
288     f2 = fCorrObj->GetVZEROC();
289     f3 = fCorrObj->GetVZEROCorr();
290   }
291   
292   (*fpcstream) <<"event"<<
293     "nTracksTPC=" << i1 <<
294     "nTracks="    << i2 <<
295     "v0A="        << f1 <<
296     "v0C="        << f2 <<
297     "v0Corr="     << f3 << "\n";
298 #endif
299
300   // -- Post output data
301   // ---------------------
302   PostData(1,fOutList);
303
304   return;
305 }      
306
307 /*
308  * ---------------------------------------------------------------------------------
309  *                            Setup Methods - private
310  * ---------------------------------------------------------------------------------
311  */
312
313 //________________________________________________________________________
314 Bool_t AliAnalysisTaskHIMultCorr::SetupEvent() {
315   // Setup Reading of event
316   
317   Bool_t aCuts[] = {kFALSE, kFALSE, kFALSE, kFALSE, kFALSE, kFALSE, kFALSE};
318
319   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
320   // -- ESD Event
321   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
322
323   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> 
324     (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
325   if (!esdH) {
326     AliError("Could not get ESDInputHandler");
327     return kFALSE;
328   } 
329
330   fESD = esdH->GetEvent();
331   if (!fESD) {
332     AliError("fESD not available");
333     return kFALSE;
334   }
335
336   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
337   // -- Check object existence
338   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
339
340   const AliESDVertex*    vtxESD    = fESD->GetPrimaryVertexTracks();
341   const AliESDVertex*    vtxESDTPC = fESD->GetPrimaryVertexTPC();  
342   const AliESDVertex*    vtxESDSPD = fESD->GetPrimaryVertexSPD();  
343   const AliMultiplicity* multESD   = fESD->GetMultiplicity();  
344
345   if ( !vtxESD ){
346     AliError("No Tracks Vertex");
347     return kFALSE;
348   }
349
350   if ( !vtxESDTPC ){
351     AliError("No TPC Vertex");
352     return kFALSE;
353   }
354
355   if ( !vtxESDSPD ){
356     AliError("No SPD Vertex");
357     return kFALSE;
358   }
359
360   if ( !multESD ){
361     AliError("No Multiplicity");
362     return kFALSE;
363   }
364
365   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
366   // -- After Physics Selection;
367   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
368
369   aCuts[0] = kFALSE;
370
371   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
372   // -- ZDC cut
373   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
374
375   if ( !fIsMC ) {
376     if (!fTriggerAnalysis->ZDCTimeTrigger(fESD)) 
377       aCuts[1] = kTRUE;
378   }
379
380   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
381   // -- Francesco Prino cut's
382   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
383
384   if (!vtxESDTPC->GetStatus()) 
385     aCuts[2] = kTRUE;
386
387   if( vtxESDTPC->GetNContributors() < (-10.+0.25*multESD->GetNumberOfITSClusters(0)) ) 
388     aCuts[3] = kTRUE;
389
390   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
391   // -- Vertex Cuts - Tracks 
392   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
393
394   //  if (!vtxESD->GetStatus())
395   //    aCuts[4] = kTRUE;
396   
397   if(vtxESDTPC->GetZv() > fMaxVertexZ || vtxESDTPC->GetZv() < (-1.*fMaxVertexZ)) 
398     aCuts[5] = kTRUE;
399
400   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
401   // -- CentralityBin
402   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
403
404   if (GetCentralityBin() < 0)
405     aCuts[6] = kTRUE;
406
407   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
408   // -- Marian's debug streamer
409   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
410 #if USE_STREAMER
411   //  TFile *inFile=(dynamic_cast<TChain*>(GetInputData(0)))->GetCurrentFile();
412   //  TString fileName(inFile->GetName()); 
413   static Int_t runNo = fESD->GetRunNumber();
414   static  UInt_t tt = fESD->GetTimeStamp();
415   static Float_t zvTPC = vtxESDTPC->GetZv();
416   static Float_t zvSPD = vtxESDSPD->GetZv();
417   static Float_t zvTra = vtxESD->GetZv();
418   static Int_t ncontTPCV = vtxESDTPC->GetNContributors();
419   static Int_t ncontSPDV = vtxESDSPD->GetNContributors();
420   static Int_t ncontTraV = vtxESD->GetNContributors();
421   static Bool_t cutPrino = aCuts[3];
422   static Bool_t cutZDCTi = aCuts[1];
423   static Int_t spd1 =  multESD->GetNumberOfITSClusters(0);
424   static Int_t spd2 =  multESD->GetNumberOfITSClusters(1);
425
426   (*fpcstream) << "event" <<
427     "run="     << runNo                << //runNo
428     "time="    << tt                   <<
429     //    "fname="   << fileName             <<
430     //    "eventNr=" << eventNumber<<
431     "zvTPC="   << zvTPC                << //zvertex TPC
432     "zvSPD="   << zvSPD                << //zvertex SPD
433     "zvTra="   << zvTra                << //zvertex Tracks
434     "ncontTPCV"<< ncontTPCV            << // N contributors to TPC vtx
435     "ncontSPDV"<< ncontSPDV            << // N contributors to SPD vtx
436     "ncontTraV"<< ncontTraV            << // N contributors to Tracks vtx
437     "cutPrino="<< cutPrino             << //francecos cut 
438     "cutZDCTi="<< cutZDCTi             << //zdctiming cut 
439     "spd1="    << spd1                 << // N clus in SPD0
440     "spd2="    << spd2                 << // N clus in SPD1
441     "centSPD=" << fCentralitySPD       <<
442     "centSPDB="<< fCentralitySPDBin    <<
443     "centV0="  << fCentralityVZERO     <<
444     "centV0B=" << fCentralityVZEROBin;
445 #endif
446   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
447   // -- Fill statistics / reject event
448   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
449  
450   Bool_t isRejected = kFALSE;
451
452   for (Int_t idx = 0; idx < 7; ++idx) {
453     if (aCuts[idx])
454       isRejected = kTRUE;
455     else
456       fHStat->Fill(idx);
457   }
458   
459   for (Int_t idx = 0; idx < 7; ++idx) {
460     if (aCuts[idx])
461       break;
462     fHStat->Fill(10+idx);
463   }
464
465   if (isRejected)
466     return kFALSE;
467   
468   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
469   // -- Fill Centrality Histograms
470   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
471
472   (static_cast<TH1F*>(fOutList->FindObject("fVZEROCentrality")))->Fill(fCentralityVZERO);
473   (static_cast<TH1F*>(fOutList->FindObject("fVZEROCentralityBin")))->Fill(fCentralityVZEROBin);
474
475   (static_cast<TH1F*>(fOutList->FindObject("fSPDCentrality")))->Fill(fCentralitySPD);
476   (static_cast<TH1F*>(fOutList->FindObject("fSPDCentralityBin")))->Fill(fCentralitySPDBin);
477
478   (static_cast<TH2F*>(fOutList->FindObject("fCorrCentrality")))->Fill(fCentralitySPD,fCentralityVZERO);
479   (static_cast<TH2F*>(fOutList->FindObject("fCorrCentralityBin")))->Fill(fCentralitySPDBin,fCentralityVZEROBin);
480
481   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
482   // -- Fill Vertex Histograms
483   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
484   
485   (static_cast<TH2F*>(fOutList->FindObject("fCorrVertexTracks")))->Fill(vtxESD->GetStatus(),vtxESD->GetNContributors());
486   (static_cast<TH2F*>(fOutList->FindObject("fCorrVertexTPC")))->Fill(vtxESDTPC->GetStatus(),vtxESDTPC->GetNContributors());
487   (static_cast<TH2F*>(fOutList->FindObject("fCrossCorrVertexTracks")))->Fill(vtxESDTPC->GetStatus(),vtxESD->GetNContributors());
488   (static_cast<TH2F*>(fOutList->FindObject("fCrossCorrVertexTPC")))->Fill(vtxESD->GetStatus(),vtxESDTPC->GetNContributors());
489   (static_cast<TH2F*>(fOutList->FindObject("fCorrVertexNCont")))->Fill(vtxESD->GetNContributors(),vtxESDTPC->GetNContributors());
490
491   (static_cast<TH1F*>(fOutList->FindObject("fDeltaVx")))->Fill(vtxESD->GetXv() - vtxESDTPC->GetXv());
492   (static_cast<TH1F*>(fOutList->FindObject("fDeltaVy")))->Fill(vtxESD->GetYv() - vtxESDTPC->GetYv());
493   (static_cast<TH1F*>(fOutList->FindObject("fDeltaVz")))->Fill(vtxESD->GetZv() - vtxESDTPC->GetZv());
494
495   (static_cast<TH2F*>(fOutList->FindObject("fVzTPCvsSPD")))->Fill(vtxESDTPC->GetZv(),vtxESDSPD->GetZv());
496   (static_cast<TH2F*>(fOutList->FindObject("fVzTPCvsTracks")))->Fill(vtxESDTPC->GetZv(),vtxESD->GetZv());
497   (static_cast<TH2F*>(fOutList->FindObject("fVzSPDvsTracks")))->Fill(vtxESDSPD->GetZv(),vtxESD->GetZv());
498
499   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --   
500
501   return kTRUE;
502 }
503
504 //________________________________________________________________________
505 Int_t AliAnalysisTaskHIMultCorr::GetCentralityBin(){
506   // Get Centrality bin
507   
508   fCentralityBin = -1;
509
510   if (fUseCentralitySel == 0)
511     return fCentralityBin;
512   
513   AliCentrality *esdCentrality = fESD->GetCentrality();
514   fCentralityVZERO  = esdCentrality->GetCentralityPercentile("V0M");  
515   fCentralitySPD    = esdCentrality->GetCentralityPercentile("CL1");
516
517   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --
518
519   fCentralityVZEROBin = esdCentrality->GetCentralityClass10("V0M")*10;
520   if (fCentralityVZEROBin == 0 && esdCentrality->GetCentralityClass5("V0M") == 1) 
521     fCentralityVZEROBin = 5;
522
523   fCentralitySPDBin = esdCentrality->GetCentralityClass10("CL1")*10;
524   if (fCentralitySPDBin == 0 && esdCentrality->GetCentralityClass5("CL1") == 1) 
525     fCentralitySPDBin = 5;
526
527   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --
528
529   if ( fUseCentralitySel == 1 )
530     fCentralityBin = fCentralityVZEROBin;
531   else if ( fUseCentralitySel == 2 )
532     fCentralityBin = fCentralitySPDBin;
533
534   if ( fCentralityBin == 100 )
535     fCentralityBin = -1;
536
537   // -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- --
538       
539   return fCentralityBin;
540 }
541
542 //________________________________________________________________________
543 void AliAnalysisTaskHIMultCorr::Terminate(Option_t *){
544   // Terminate
545 #if USE_STREAMER
546   if (fpcstream) fpcstream->GetFile()->Write();
547 #endif
548
549 }