]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/physics/AliHLTMultiplicityCorrelations.cxx
Reading friends in analysis framework inside HLT
[u/mrichter/AliRoot.git] / HLT / global / physics / AliHLTMultiplicityCorrelations.cxx
1 //-*- Mode: C++ -*-
2 // $Id: AliHLTMultiplicityCorrelations.cxx  $
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Jochen Thaeder <jochen@thaeder.de>                    *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTMultiplicityCorrelations.cxx
20     @author Jochen Thaeder <jochen@thaeder.de>
21     @brief  Correlation plots for multiplicity studies
22 */
23
24 // see header file for class documentation
25 // or
26 // refer to README to build package
27 // or
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
29
30 #include "TH1F.h"
31 #include "TH2F.h"
32 #include "TList.h"
33
34 #include "AliHLTMultiplicityCorrelations.h"
35
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTMultiplicityCorrelations)
38
39 #if __GNUC__>= 3
40    using namespace std;
41 #endif
42
43 /*
44  * ---------------------------------------------------------------------------------
45  *                            Constructor / Destructor
46  * ---------------------------------------------------------------------------------
47  */
48
49 //##################################################################################
50 AliHLTMultiplicityCorrelations::AliHLTMultiplicityCorrelations() :
51   TNamed(),
52   fHistList(NULL),
53   fESDEvent(NULL),
54   fESDZDC(NULL),
55   fESDVZERO(NULL),
56   fESDTrackCuts(NULL),
57   fCentHistV0Mpercentile(NULL),
58   fProcessTPC(kTRUE), fProcessSPD(kTRUE),
59   fProcessVZERO(kTRUE), fProcessZDC(kTRUE), fProcessCentrality(kTRUE),
60   fEsdTracks(0), fEsdTracksA(0),
61   fTpcTracks(0), fTpcTracksA(0),
62   fTpcTracksRef(0),
63   fVzeroMult(0.), fVzeroTriggerMult(0.),
64   fSpdNClusters(0), fSpdNClustersInner(0), fSpdNClustersOuter(0),
65   fVzeroBinning(1500), fVzeroBinningMin(0.), fVzeroBinningMax(30000.),
66   fTpcBinning(800),fTpcBinningMin(0.),fTpcBinningMax(8000.),
67   fZdcBinning(280),fZdcBinningMin(0.),fZdcBinningMax(140.),
68   fZemBinning(100),fZemBinningMin(0.),fZemBinningMax(5.),
69   fSpdBinning(750),fSpdBinningMin(0.),fSpdBinningMax(15000.) {
70   // see header file for class documentation
71   // or
72   // refer to README to build package
73   // or
74   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
75   
76 }
77 //##################################################################################
78 AliHLTMultiplicityCorrelations::AliHLTMultiplicityCorrelations(Char_t* name, Char_t* title) : 
79   TNamed(name,title), 
80   fHistList(NULL),
81   fESDEvent(NULL),
82   fESDZDC(NULL),
83   fESDVZERO(NULL),
84   fESDTrackCuts(NULL),
85   fCentHistV0Mpercentile(NULL),
86   fProcessTPC(kTRUE), fProcessSPD(kTRUE),
87   fProcessVZERO(kTRUE), fProcessZDC(kTRUE), fProcessCentrality(kTRUE),
88   fEsdTracks(0), fEsdTracksA(0),
89   fTpcTracks(0), fTpcTracksA(0),
90   fTpcTracksRef(0),
91   fVzeroMult(0.), fVzeroTriggerMult(0.),
92   fSpdNClusters(0), fSpdNClustersInner(0), fSpdNClustersOuter(0),
93   fVzeroBinning(1500), fVzeroBinningMin(0.), fVzeroBinningMax(30000.),
94   fTpcBinning(800),fTpcBinningMin(0.),fTpcBinningMax(8000.),
95   fZdcBinning(280),fZdcBinningMin(0.),fZdcBinningMax(140.),
96   fZemBinning(100),fZemBinningMin(0.),fZemBinningMax(5.),
97   fSpdBinning(750),fSpdBinningMin(0.),fSpdBinningMax(15000.) {
98   // see header file for class documentation
99   // or
100   // refer to README to build package
101   // or
102   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
103   
104 }
105
106 //##################################################################################
107 AliHLTMultiplicityCorrelations::~AliHLTMultiplicityCorrelations() {
108   // see header file for class documentation
109
110   if ( fHistList ) {
111     fHistList->Clear();
112     delete fHistList;
113   }
114   fHistList = NULL;
115 }
116
117 /*
118  * ---------------------------------------------------------------------------------
119  *                                   Initialize / Reset
120  * ---------------------------------------------------------------------------------
121  */
122
123 //##################################################################################
124 Int_t AliHLTMultiplicityCorrelations::Initialize() {
125   // see header file for class documentation  
126
127   Int_t iResult = 0;
128
129   fHistList = new TList();
130   fHistList->SetOwner(kTRUE);
131   fHistList->SetName("MultiplicityCorrelations");
132   iResult = SetupHistograms();
133   
134   if (fProcessZDC)   { HLTInfo ("Processing of ZDC enabled"); }
135   if (fProcessTPC)   { HLTInfo ("Processing of TPC enabled"); }
136   if (fProcessSPD)   { HLTInfo ("Processing of SPD enabled"); }
137   if (fProcessVZERO) { HLTInfo ("Processing of VZERO enabled"); }
138   if (fProcessCentrality) { HLTInfo ("Processing of Centrality enabled"); }
139
140   return iResult;
141 }
142
143 //##################################################################################
144 Int_t AliHLTMultiplicityCorrelations::Initialize( const Char_t* listName) {
145   // see header file for class documentation  
146
147   Int_t iResult = 0;
148
149   fHistList = new TList();
150   fHistList->SetOwner(kTRUE);
151   fHistList->SetName(Form("MultiplicityCorrelations_%s",listName));
152   iResult = SetupHistograms();
153   
154   if (fProcessZDC)   { HLTInfo ("Processing of ZDC enabled"); }
155   if (fProcessTPC)   { HLTInfo ("Processing of TPC enabled"); }
156   if (fProcessSPD)   { HLTInfo ("Processing of SPD enabled"); }
157   if (fProcessVZERO) { HLTInfo ("Processing of VZERO enabled"); }
158   if (fProcessCentrality) { HLTInfo ("Processing of Centrality enabled"); }
159
160   return iResult;
161 }
162
163 /*
164  * ---------------------------------------------------------------------------------
165  *                             Output - public
166  * ---------------------------------------------------------------------------------
167  */
168
169 //##################################################################################
170 Int_t AliHLTMultiplicityCorrelations::ProcessEvent( AliESDEvent *esd, AliESDVZERO* esdVZERO, Int_t nSpdClusters) {
171   // see header file for class documentation  
172
173   Int_t iResult = 0;
174
175   if ( ! AddESDEvent(esd) ) {
176     HLTWarning("No ESD event.");
177     return -1;
178   }
179   
180   if ( esdVZERO )
181     fESDVZERO = esdVZERO;
182
183   // -- TPC .. To be done before the others
184   if (fESDEvent->GetNumberOfTracks() > 0 && fProcessTPC)
185     iResult = ProcessTPC();
186   
187   // -- SPD
188   fSpdNClusters = nSpdClusters;
189   if (fProcessSPD)
190     iResult = ProcessSPD();
191     
192   // -- VZERO
193   if (fESDVZERO && fProcessVZERO)
194     iResult = ProcessVZERO();
195
196   // -- ZDC and Correlations
197   if (fESDZDC && fProcessZDC)
198     iResult = ProcessZDC();
199  
200   // -- Centrality
201   if (fCentHistV0Mpercentile && fProcessCentrality && fVzeroMult > 0.)
202     iResult = ProcessCentrality();
203
204   return iResult;
205 }
206
207 ////////////////////////////////////////////////////////////////////////////////////
208 ////////////////////////////////////////////////////////////////////////////////////
209 ///                                                                              ///
210 //////                             PRIVATE                                    //////
211 ///                                                                              ///
212 ////////////////////////////////////////////////////////////////////////////////////
213 ////////////////////////////////////////////////////////////////////////////////////
214
215 /*
216  * ---------------------------------------------------------------------------------
217  *                          Setup / Initialize - private
218  * ---------------------------------------------------------------------------------
219  */
220
221 //##################################################################################
222 Bool_t AliHLTMultiplicityCorrelations::AddESDEvent( AliESDEvent* esd ) {
223   // see header file for class documentation  
224
225   fESDEvent = esd;
226   
227   // -- Check for ESD
228   if ( !fESDEvent ) {
229     HLTWarning("No ESD event present.");
230     return kFALSE;
231   }
232   
233   // -- Check for PrimaryVertex
234   if ( !esd->GetPrimaryVertexTracks() ){
235     HLTError("No Vertex present.");
236     return kFALSE;
237   }
238   
239   fESDZDC = esd->GetESDZDC();
240   if ( !fESDZDC ) {
241     HLTInfo("No ZDC information !");
242   }
243   
244   fESDVZERO = esd->GetVZEROData();
245   if ( !fESDVZERO ) {
246     HLTInfo("No VZERO information !");
247   }
248   
249   return kTRUE;
250 }
251
252 //##################################################################################
253 Int_t AliHLTMultiplicityCorrelations::SetupHistograms() {
254   // see header file for class documentation  
255
256   Int_t iResult = 0;
257   
258   if (fProcessVZERO) iResult = SetupVZERO();
259   if (fProcessZDC)   iResult = SetupZDC();
260   if (fProcessTPC)   iResult = SetupTPC();
261   if (fProcessSPD)   iResult = SetupSPD();
262   if (fProcessCentrality)   iResult = SetupCentrality();
263
264   iResult = SetupCorrelations();
265
266   return iResult;
267 }
268
269 //##################################################################################
270 Int_t AliHLTMultiplicityCorrelations::SetupVZERO() {
271   // see header file for class documentation  
272
273   // VzeroMult
274   fHistList->Add(new TH1F("fVzeroMult",   "Multiplicity^{VZERO};Multiplicity^{VZERO};N_{Events}",   
275                           fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax));
276   fHistList->Add(new TH2F("fVzeroMultAC", "Multiplicity^{VZERO} A vs C;Multiplicity^{VZERO}A ;Multiplicity^{VZERO} C", 
277                           fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax)); 
278
279   // VzeroTriggerMult 
280   fHistList->Add(new TH1F("fVzeroTriggerMult",   "Trigger Multiplicity^{VZERO};Trigger Multiplicity^{VZERO};N_{Events}",   
281                           fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax)); 
282   fHistList->Add(new TH2F("fVzeroTriggerMultAC", "Trigger Multiplicity^{VZERO} A vs C;Trigger Multiplicity^{VZERO}A ;Trigger Multiplicity^{VZERO} C", 
283                           fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax)); 
284
285   return 0;
286 }
287
288 //##################################################################################
289 Int_t AliHLTMultiplicityCorrelations::SetupZDC() {
290   // see header file for class documentation  
291   
292   // E_{ZDC}
293   fHistList->Add(new TH1F("fZdcEzdc",  "E_{ZDC};E_{ZDC} (TeV);N_{Events}",   
294                           fZdcBinning,fZdcBinningMin,fZdcBinningMax));
295   fHistList->Add(new TH2F("fZdcEzdcAEzdcC", "E_{ZDC} A vs E_{ZDC} C;E_{ZDC} A (TeV); E_{ZDC} C (TeV)",   
296                           fZdcBinning,fZdcBinningMin,fZdcBinningMax, fZdcBinning,fZdcBinningMin,fZdcBinningMax));
297
298   // E_{ZEM} vs E_{ZDC} 
299   fHistList->Add(new TH2F("fZdcEzemEzdc", "E_{ZEM} vs E_{ZDC};E_{ZEM} (TeV); E_{ZDC} (TeV)",   
300                           fZemBinning,fZemBinningMin,fZemBinningMax, fZdcBinning,fZdcBinningMin,fZdcBinningMax));
301
302   return 0;
303 }
304
305 //##################################################################################
306 Int_t AliHLTMultiplicityCorrelations::SetupTPC() {
307   // see header file for class documentation  
308
309   // Multiplicty
310   fHistList->Add(new TH1F("fTpcNch0", "TPC N_{ch} esdTracks; TPC N_{ch};N_{Events}", 
311                           fTpcBinning,fTpcBinningMin,fTpcBinningMax));
312   fHistList->Add(new TH1F("fTpcNch1", "TPC N_{ch} accepted esdTracks; TPC N_{ch};N_{Events}", 
313                           fTpcBinning,fTpcBinningMin,fTpcBinningMax));
314   fHistList->Add(new TH1F("fTpcNch2", "TPC N_{ch} tpcTracks; TPC N_{ch};N_{Events}", 
315                           fTpcBinning,fTpcBinningMin,fTpcBinningMax));
316   fHistList->Add(new TH1F("fTpcNch3", "TPC N_{ch} accepted tpcTracks; TPC N_{ch};N_{Events}", 
317                           fTpcBinning,fTpcBinningMin,fTpcBinningMax));
318   fHistList->Add(new TH1F("fTpcNch4", "TPC N_{ch} reference Tracks; TPC N_{ch};N_{Events}", 
319                           fTpcBinning,fTpcBinningMin,fTpcBinningMax));
320
321   return 0;
322 }
323
324 //##################################################################################
325 Int_t AliHLTMultiplicityCorrelations::SetupCorrelations() {
326   // see header file for class documentation  
327
328   // ----------------------------------------------------
329   //   ZDC vs TPC
330   // ----------------------------------------------------
331   if (fProcessTPC && fProcessZDC) {
332
333     // E_{ZDC} vs N_{ch}
334     fHistList->Add(new TH2F("fCorrEzdcTpc",    "E_{ZDC} vs N_{ch}^{TPC};E_{ZDC} (TeV);N_{ch}^{TPC}",   
335                             fZdcBinning,fZdcBinningMin,fZdcBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));
336
337     fHistList->Add(new TH2F("fCorrEzdcRefTpc", "E_{ZDC} vs N_{ch}^{Ref TPC};E_{ZDC} (TeV);N_{ch}^{Ref TPC}",   
338                             fZdcBinning,fZdcBinningMin,fZdcBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));
339   }
340   // ----------------------------------------------------
341   //   ZDC vs VZERO
342   // ----------------------------------------------------
343   if (fProcessZDC && fProcessVZERO) {
344
345     // E_{ZDC} vs Multiplicty VZERO
346     fHistList->Add(new TH2F("fCorrEzdcVzero", 
347                             "E_{ZDC} vs Multiplicity^{VZERO};E_{ZDC} (TeV);Multiplicity^{VZERO}",  
348                             fZdcBinning,fZdcBinningMin,fZdcBinningMax, fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax));
349     fHistList->Add(new TH2F("fCorrEzdcTriggerVzero", 
350                             "E_{ZDC} vs Trigger Multiplicity^{VZERO};E_{ZDC} (TeV);Trigger Multiplicity^{VZERO}",  
351                             fZdcBinning,fZdcBinningMin,fZdcBinningMax, fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax));
352   }
353   // ----------------------------------------------------
354   //   VZERO vs TPC
355   // ----------------------------------------------------
356   if ( fProcessTPC && fProcessVZERO ) {
357
358     fHistList->Add(new TH2F("fCorrVzeroTpc", 
359                             "Multiplicity^{VZERO} vs N_{ch}^{TPC};Multiplicity^{VZERO};N_{ch}^{TPC}", 
360                             fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));
361     fHistList->Add(new TH2F("fCorrVzeroRefTpc", 
362                             "Multiplicity^{VZERO} vs N_{ch}^{Ref TPC};Multiplicity^{VZERO};N_{ch}^{Ref TPC}", 
363                             fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));
364   
365     fHistList->Add(new TH2F("fCorrTriggerVzeroTpc", 
366                             "Trigger Multiplicity^{VZERO} vs N_{ch}^{TPC};Multiplicity^{VZERO};N_{ch}^{TPC}", 
367                             fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));
368     fHistList->Add(new TH2F("fCorrTriggerVzeroRefTpc", 
369                             "Trigger Multiplicity^{VZERO} vs N_{ch}^{Ref TPC};Multiplicity^{VZERO};N_{ch}^{Ref TPC}", 
370                             fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));
371   }
372   // ----------------------------------------------------
373   //   SPD vs TPC
374   // ----------------------------------------------------
375   if ( fProcessSPD && fProcessTPC ) {
376     fHistList->Add(new TH2F("fCorrSpdOuterTpc", "N_{clusters}^{SPD}_{Outer} vs N_{ch}^{TPC};N_{clusters}^{SPD};N_{ch}^{TPC}",   
377                             fSpdBinning,fSpdBinningMin,fSpdBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));  
378     fHistList->Add(new TH2F("fCorrSpdOuterRefTpc", "N_{clusters}^{SPD}_{Outer} vs N_{ch}^{Ref TPC};N_{clusters}^{SPD};N_{ch}^{Ref TPC}",   
379                             fSpdBinning,fSpdBinningMin,fSpdBinningMax, fTpcBinning,fTpcBinningMin,fTpcBinningMax));  
380     
381
382   }
383   // ----------------------------------------------------
384   //   SPD vs VZERO
385   // ----------------------------------------------------
386   if ( fProcessSPD && fProcessVZERO ) {
387     fHistList->Add(new TH2F("fCorrVzeroSpdOuter", 
388                             "Multiplicity^{VZERO} vs N_{ch}^{SPD}_{Outer};Multiplicity^{VZERO};N^{SPD}", 
389                             fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fSpdBinning,fSpdBinningMin,fSpdBinningMax));
390     
391     fHistList->Add(new TH2F("fCorrVzeroTriggerSpdOuter", 
392                             "Trigger Multiplicity^{VZERO} vs N_{ch}^{SPD}_{Outer};Trigger Multiplicity^{VZERO};N^{SPD}", 
393                             fVzeroBinning,fVzeroBinningMin,fVzeroBinningMax, fSpdBinning,fSpdBinningMin,fSpdBinningMax));
394     
395   }
396   // ----------------------------------------------------
397   //   SPD vs ZDC
398   // ----------------------------------------------------
399   if ( fProcessSPD && fProcessZDC) {
400     fHistList->Add(new TH2F("fCorrEzdcSpdOuter", "E_{ZDC} vs N_{ch}^{SPD}_{Outer};E_{ZDC} (TeV);N^{SPD}",   
401                             fZdcBinning,fZdcBinningMin,fZdcBinningMax, fSpdBinning,fSpdBinningMin,fSpdBinningMax));
402   }
403
404   return 0;
405 }
406
407 //##################################################################################
408 Int_t AliHLTMultiplicityCorrelations::SetupSPD() {
409   // see header file for class documentation  
410   
411   fHistList->Add(new TH1F("fSpdNClusters", "Multplicity_{SPD};Multplicity_{SPD};N_{Events}",   
412                           fSpdBinning,fSpdBinningMin,fSpdBinningMax));
413
414   fHistList->Add(new TH1F("fSpdNClustersInner", "Multplicity_{SPD} Layer 0;Multplicity_{SPD};N_{Events}",   
415                           fSpdBinning,fSpdBinningMin,fSpdBinningMax));
416
417   fHistList->Add(new TH1F("fSpdNClustersOuter", "Multplicity_{SPD} Layer 1;Multplicity_{SPD};N_{Events}",   
418                           fSpdBinning,fSpdBinningMin,fSpdBinningMax));
419
420   return 0;
421 }
422 //##################################################################################
423 Int_t AliHLTMultiplicityCorrelations::SetupCentrality() {
424   // see header file for class documentation  
425
426   
427   if (!fCentHistV0Mpercentile) {
428     fProcessCentrality = kFALSE;
429     return -1;  
430   }
431
432   // Vzero Centrality
433   fHistList->Add(new TH1F("fVzeroCentV0M", "Centrality V0M;Centrality", 100, 0, 100));
434
435   return 0;
436 }
437
438 /*
439  * ---------------------------------------------------------------------------------
440  *                               Process - private
441  * ---------------------------------------------------------------------------------
442  */
443
444 //##################################################################################
445 Int_t AliHLTMultiplicityCorrelations::ProcessTPC() {
446   // see header file for class documentation  
447
448   Int_t iResult = 0 ;
449
450   fEsdTracks = 0;  
451   fEsdTracksA = 0;  
452   fTpcTracks = 0;  
453   fTpcTracksA = 0;  
454   fTpcTracksRef = 0;  
455   
456   // -----------------------------------------------------------------
457   for (Int_t idx = 0; idx < fESDEvent->GetNumberOfTracks(); idx++) {
458     AliESDtrack* esdTrack = fESDEvent->GetTrack(idx);
459     if (!esdTrack) continue;
460     ++fEsdTracks;
461     
462     if (!fESDTrackCuts->AcceptTrack(esdTrack)) continue;
463     ++fEsdTracksA;
464   }
465   // -----------------------------------------------------------------
466   for (Int_t idx = 0; idx < fESDEvent->GetNumberOfTracks(); idx++) {
467     AliESDtrack* esdTrack = fESDEvent->GetTrack(idx);
468     if (!esdTrack) continue;
469
470     // -- TPC only
471     const AliExternalTrackParam *tpcTrack = esdTrack->GetTPCInnerParam();
472     if (!tpcTrack) continue;
473     ++fTpcTracks;
474
475
476     if (!fESDTrackCuts->AcceptTrack(esdTrack)) continue;
477     ++fTpcTracksA;
478   }
479   // -----------------------------------------------------------------
480
481   fTpcTracksRef = fESDTrackCuts->GetReferenceMultiplicity(fESDEvent,kTRUE);
482
483   (static_cast<TH1F*>(fHistList->FindObject("fTpcNch0")))->Fill(fEsdTracks);
484   (static_cast<TH1F*>(fHistList->FindObject("fTpcNch1")))->Fill(fEsdTracksA);
485   (static_cast<TH1F*>(fHistList->FindObject("fTpcNch2")))->Fill(fTpcTracks);
486   (static_cast<TH1F*>(fHistList->FindObject("fTpcNch3")))->Fill(fTpcTracksA);
487   (static_cast<TH1F*>(fHistList->FindObject("fTpcNch4")))->Fill(fTpcTracksRef);
488
489   return iResult;
490 }
491
492 //##################################################################################
493 Int_t AliHLTMultiplicityCorrelations::ProcessVZERO() {
494   // see header file for class documentation  
495
496   Int_t iResult = 0 ;
497
498   Float_t vzeroMultA = fESDVZERO->GetMTotV0A();
499   Float_t vzeroMultC = fESDVZERO->GetMTotV0C();
500   fVzeroMult = vzeroMultA + vzeroMultC;
501
502   (static_cast<TH1F*>(fHistList->FindObject("fVzeroMult")))->Fill(fVzeroMult);
503   (static_cast<TH1F*>(fHistList->FindObject("fVzeroMultAC")))->Fill(vzeroMultA,vzeroMultC);
504
505   Float_t vzeroTriggerMultA = Float_t(fESDVZERO->GetTriggerChargeA());
506   Float_t vzeroTriggerMultC = Float_t(fESDVZERO->GetTriggerChargeC());
507   fVzeroTriggerMult  = vzeroTriggerMultA + vzeroTriggerMultC;
508
509   (static_cast<TH1F*>(fHistList->FindObject("fVzeroTriggerMult")))->Fill(fVzeroTriggerMult);
510   (static_cast<TH1F*>(fHistList->FindObject("fVzeroTriggerMultAC")))->Fill(vzeroTriggerMultA,vzeroTriggerMultC);
511
512   // -- VZERO - TPC correlations
513   if (fESDEvent->GetNumberOfTracks() > 0 && fProcessTPC) {
514     (static_cast<TH2F*>(fHistList->FindObject("fCorrVzeroTpc")))->Fill(fVzeroMult, fTpcTracksA);
515     (static_cast<TH2F*>(fHistList->FindObject("fCorrVzeroRefTpc")))->Fill(fVzeroMult, fTpcTracksRef);
516
517     (static_cast<TH2F*>(fHistList->FindObject("fCorrTriggerVzeroTpc")))->Fill(fVzeroTriggerMult, fTpcTracksA);
518     (static_cast<TH2F*>(fHistList->FindObject("fCorrTriggerVzeroRefTpc")))->Fill(fVzeroTriggerMult, fTpcTracksRef);
519   }
520
521   // -- VZERO - SPD correlations
522   if (fESDEvent->GetNumberOfTracks() > 0 && fProcessSPD) {
523     (static_cast<TH2F*>(fHistList->FindObject("fCorrVzeroSpdOuter")))->Fill(fVzeroMult, fSpdNClustersOuter);
524     (static_cast<TH2F*>(fHistList->FindObject("fCorrVzeroTriggerSpdOuter")))->Fill(fVzeroTriggerMult, fSpdNClustersOuter);
525   }
526   
527   return iResult;
528 }
529
530 //##################################################################################
531 Int_t AliHLTMultiplicityCorrelations::ProcessZDC() {
532   // see header file for class documentation  
533
534   Int_t iResult = 0 ;
535
536   Double_t zdcEznA = fESDZDC->GetZDCN2Energy() / 1000.;
537   Double_t zdcEznC = fESDZDC->GetZDCN1Energy() / 1000.;
538   
539   Double_t zdcEzpA = fESDZDC->GetZDCP2Energy() / 1000.;
540   Double_t zdcEzpC = fESDZDC->GetZDCP1Energy() / 1000.;
541   
542   Double_t zdcEA = (zdcEznA + zdcEzpA);
543   Double_t zdcEC = (zdcEznC + zdcEzpC);
544   Double_t zdcE = zdcEA + zdcEC;
545   
546   (static_cast<TH1F*>(fHistList->FindObject("fZdcEzdc")))->Fill(zdcE);
547   (static_cast<TH1F*>(fHistList->FindObject("fZdcEzdcAEzdcC")))->Fill(zdcEA, zdcEC);
548
549   Double_t zdcEzemA = fESDZDC->GetZDCEMEnergy(1) / 1000.;
550   Double_t zdcEzemC = fESDZDC->GetZDCEMEnergy(0) / 1000.;
551   Double_t zdcEzem  = zdcEzemA + zdcEzemC;
552   
553   (static_cast<TH2F*>(fHistList->FindObject("fZdcEzemEzdc")))->Fill(zdcEzem, zdcE);
554  
555   // -- ZDC - TPC correlations
556   if (fESDEvent->GetNumberOfTracks() > 0 && fProcessTPC) {
557     (static_cast<TH2F*>(fHistList->FindObject("fCorrEzdcTpc")))->Fill(zdcE, fTpcTracksA);
558     (static_cast<TH2F*>(fHistList->FindObject("fCorrEzdcRefTpc")))->Fill(zdcE, fTpcTracksRef);
559   }
560   
561   // -- ZDC - SPD correlations
562   if (fProcessSPD) {
563     (static_cast<TH2F*>(fHistList->FindObject("fCorrEzdcSpdOuter")))->Fill(zdcE, fSpdNClustersOuter);
564   }
565
566   // -- VZERO - ZDC correlations
567   if (fESDVZERO && fProcessVZERO) {
568     (static_cast<TH2F*>(fHistList->FindObject("fCorrEzdcVzero")))->Fill(zdcE, fVzeroMult);
569     (static_cast<TH2F*>(fHistList->FindObject("fCorrEzdcTriggerVzero")))->Fill(zdcE, fVzeroTriggerMult);
570   }
571
572   return iResult;
573 }
574
575 //##################################################################################
576 Int_t AliHLTMultiplicityCorrelations::ProcessSPD() {
577   // see header file for class documentation
578   
579   (static_cast<TH1F*>(fHistList->FindObject("fSpdNClusters")))->Fill(fSpdNClusters);
580   (static_cast<TH1F*>(fHistList->FindObject("fSpdNClustersInner")))->Fill(fSpdNClustersInner);
581   (static_cast<TH1F*>(fHistList->FindObject("fSpdNClustersOuter")))->Fill(fSpdNClustersOuter);
582
583   
584   // -- SPD vs TPC correlations
585   if (fProcessTPC) {
586     (static_cast<TH2F*>(fHistList->FindObject("fCorrSpdOuterTpc")))->Fill(fSpdNClustersOuter,fTpcTracksA);
587     (static_cast<TH2F*>(fHistList->FindObject("fCorrSpdOuterRefTpc")))->Fill(fSpdNClustersOuter,fTpcTracksRef);
588   }
589
590   return 0;
591 }
592
593 //##################################################################################
594 Int_t AliHLTMultiplicityCorrelations::ProcessCentrality() {
595   // see header file for class documentation
596
597   Double_t centV0M = fCentHistV0Mpercentile->GetBinContent(fCentHistV0Mpercentile->FindBin((fVzeroMult)));
598   (static_cast<TH1F*>(fHistList->FindObject("fVzeroCentV0M")))->Fill(centV0M);
599
600   return 0;
601 }