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. *
7 //* Primary Authors: Jochen Thaeder <jochen@thaeder.de> *
8 //* for The ALICE HLT Project. *
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 //**************************************************************************
19 /** @file AliHLTMultiplicityCorrelations.cxx
20 @author Jochen Thaeder <jochen@thaeder.de>
21 @brief Correlation plots for multiplicity studies
24 // see header file for class documentation
26 // refer to README to build package
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
34 #include "AliHLTMultiplicityCorrelations.h"
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTMultiplicityCorrelations)
44 * ---------------------------------------------------------------------------------
45 * Constructor / Destructor
46 * ---------------------------------------------------------------------------------
49 //##################################################################################
50 AliHLTMultiplicityCorrelations::AliHLTMultiplicityCorrelations() :
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),
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
72 // refer to README to build package
74 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
77 //##################################################################################
78 AliHLTMultiplicityCorrelations::AliHLTMultiplicityCorrelations(Char_t* name, Char_t* title) :
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),
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
100 // refer to README to build package
102 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
106 //##################################################################################
107 AliHLTMultiplicityCorrelations::~AliHLTMultiplicityCorrelations() {
108 // see header file for class documentation
118 * ---------------------------------------------------------------------------------
120 * ---------------------------------------------------------------------------------
123 //##################################################################################
124 Int_t AliHLTMultiplicityCorrelations::Initialize() {
125 // see header file for class documentation
129 fHistList = new TList();
130 fHistList->SetOwner(kTRUE);
131 fHistList->SetName("MultiplicityCorrelations");
132 iResult = SetupHistograms();
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"); }
143 //##################################################################################
144 Int_t AliHLTMultiplicityCorrelations::Initialize( const Char_t* listName) {
145 // see header file for class documentation
149 fHistList = new TList();
150 fHistList->SetOwner(kTRUE);
151 fHistList->SetName(Form("MultiplicityCorrelations_%s",listName));
152 iResult = SetupHistograms();
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"); }
164 * ---------------------------------------------------------------------------------
166 * ---------------------------------------------------------------------------------
169 //##################################################################################
170 Int_t AliHLTMultiplicityCorrelations::ProcessEvent( AliESDEvent *esd, AliESDVZERO* esdVZERO, Int_t nSpdClusters) {
171 // see header file for class documentation
175 if ( ! AddESDEvent(esd) ) {
176 HLTWarning("No ESD event.");
181 fESDVZERO = esdVZERO;
183 // -- TPC .. To be done before the others
184 if (fESDEvent->GetNumberOfTracks() > 0 && fProcessTPC)
185 iResult = ProcessTPC();
188 fSpdNClusters = nSpdClusters;
190 iResult = ProcessSPD();
193 if (fESDVZERO && fProcessVZERO)
194 iResult = ProcessVZERO();
196 // -- ZDC and Correlations
197 if (fESDZDC && fProcessZDC)
198 iResult = ProcessZDC();
201 if (fCentHistV0Mpercentile && fProcessCentrality && fVzeroMult > 0.)
202 iResult = ProcessCentrality();
207 ////////////////////////////////////////////////////////////////////////////////////
208 ////////////////////////////////////////////////////////////////////////////////////
210 ////// PRIVATE //////
212 ////////////////////////////////////////////////////////////////////////////////////
213 ////////////////////////////////////////////////////////////////////////////////////
216 * ---------------------------------------------------------------------------------
217 * Setup / Initialize - private
218 * ---------------------------------------------------------------------------------
221 //##################################################################################
222 Bool_t AliHLTMultiplicityCorrelations::AddESDEvent( AliESDEvent* esd ) {
223 // see header file for class documentation
229 HLTWarning("No ESD event present.");
233 // -- Check for PrimaryVertex
234 if ( !esd->GetPrimaryVertexTracks() ){
235 HLTError("No Vertex present.");
239 fESDZDC = esd->GetESDZDC();
241 HLTInfo("No ZDC information !");
244 fESDVZERO = esd->GetVZEROData();
246 HLTInfo("No VZERO information !");
252 //##################################################################################
253 Int_t AliHLTMultiplicityCorrelations::SetupHistograms() {
254 // see header file for class documentation
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();
264 iResult = SetupCorrelations();
269 //##################################################################################
270 Int_t AliHLTMultiplicityCorrelations::SetupVZERO() {
271 // see header file for class documentation
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));
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));
288 //##################################################################################
289 Int_t AliHLTMultiplicityCorrelations::SetupZDC() {
290 // see header file for class documentation
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));
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));
305 //##################################################################################
306 Int_t AliHLTMultiplicityCorrelations::SetupTPC() {
307 // see header file for class documentation
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));
324 //##################################################################################
325 Int_t AliHLTMultiplicityCorrelations::SetupCorrelations() {
326 // see header file for class documentation
328 // ----------------------------------------------------
330 // ----------------------------------------------------
331 if (fProcessTPC && fProcessZDC) {
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));
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));
340 // ----------------------------------------------------
342 // ----------------------------------------------------
343 if (fProcessZDC && fProcessVZERO) {
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));
353 // ----------------------------------------------------
355 // ----------------------------------------------------
356 if ( fProcessTPC && fProcessVZERO ) {
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));
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));
372 // ----------------------------------------------------
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));
383 // ----------------------------------------------------
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));
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));
396 // ----------------------------------------------------
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));
407 //##################################################################################
408 Int_t AliHLTMultiplicityCorrelations::SetupSPD() {
409 // see header file for class documentation
411 fHistList->Add(new TH1F("fSpdNClusters", "Multplicity_{SPD};Multplicity_{SPD};N_{Events}",
412 fSpdBinning,fSpdBinningMin,fSpdBinningMax));
414 fHistList->Add(new TH1F("fSpdNClustersInner", "Multplicity_{SPD} Layer 0;Multplicity_{SPD};N_{Events}",
415 fSpdBinning,fSpdBinningMin,fSpdBinningMax));
417 fHistList->Add(new TH1F("fSpdNClustersOuter", "Multplicity_{SPD} Layer 1;Multplicity_{SPD};N_{Events}",
418 fSpdBinning,fSpdBinningMin,fSpdBinningMax));
422 //##################################################################################
423 Int_t AliHLTMultiplicityCorrelations::SetupCentrality() {
424 // see header file for class documentation
427 if (!fCentHistV0Mpercentile) {
428 fProcessCentrality = kFALSE;
433 fHistList->Add(new TH1F("fVzeroCentV0M", "Centrality V0M;Centrality", 100, 0, 100));
439 * ---------------------------------------------------------------------------------
441 * ---------------------------------------------------------------------------------
444 //##################################################################################
445 Int_t AliHLTMultiplicityCorrelations::ProcessTPC() {
446 // see header file for class documentation
456 // -----------------------------------------------------------------
457 for (Int_t idx = 0; idx < fESDEvent->GetNumberOfTracks(); idx++) {
458 AliESDtrack* esdTrack = fESDEvent->GetTrack(idx);
459 if (!esdTrack) continue;
462 if (!fESDTrackCuts->AcceptTrack(esdTrack)) continue;
465 // -----------------------------------------------------------------
466 for (Int_t idx = 0; idx < fESDEvent->GetNumberOfTracks(); idx++) {
467 AliESDtrack* esdTrack = fESDEvent->GetTrack(idx);
468 if (!esdTrack) continue;
471 const AliExternalTrackParam *tpcTrack = esdTrack->GetTPCInnerParam();
472 if (!tpcTrack) continue;
476 if (!fESDTrackCuts->AcceptTrack(esdTrack)) continue;
479 // -----------------------------------------------------------------
481 fTpcTracksRef = fESDTrackCuts->GetReferenceMultiplicity(fESDEvent,kTRUE);
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);
492 //##################################################################################
493 Int_t AliHLTMultiplicityCorrelations::ProcessVZERO() {
494 // see header file for class documentation
498 Float_t vzeroMultA = fESDVZERO->GetMTotV0A();
499 Float_t vzeroMultC = fESDVZERO->GetMTotV0C();
500 fVzeroMult = vzeroMultA + vzeroMultC;
502 (static_cast<TH1F*>(fHistList->FindObject("fVzeroMult")))->Fill(fVzeroMult);
503 (static_cast<TH1F*>(fHistList->FindObject("fVzeroMultAC")))->Fill(vzeroMultA,vzeroMultC);
505 Float_t vzeroTriggerMultA = Float_t(fESDVZERO->GetTriggerChargeA());
506 Float_t vzeroTriggerMultC = Float_t(fESDVZERO->GetTriggerChargeC());
507 fVzeroTriggerMult = vzeroTriggerMultA + vzeroTriggerMultC;
509 (static_cast<TH1F*>(fHistList->FindObject("fVzeroTriggerMult")))->Fill(fVzeroTriggerMult);
510 (static_cast<TH1F*>(fHistList->FindObject("fVzeroTriggerMultAC")))->Fill(vzeroTriggerMultA,vzeroTriggerMultC);
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);
517 (static_cast<TH2F*>(fHistList->FindObject("fCorrTriggerVzeroTpc")))->Fill(fVzeroTriggerMult, fTpcTracksA);
518 (static_cast<TH2F*>(fHistList->FindObject("fCorrTriggerVzeroRefTpc")))->Fill(fVzeroTriggerMult, fTpcTracksRef);
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);
530 //##################################################################################
531 Int_t AliHLTMultiplicityCorrelations::ProcessZDC() {
532 // see header file for class documentation
536 Double_t zdcEznA = fESDZDC->GetZDCN2Energy() / 1000.;
537 Double_t zdcEznC = fESDZDC->GetZDCN1Energy() / 1000.;
539 Double_t zdcEzpA = fESDZDC->GetZDCP2Energy() / 1000.;
540 Double_t zdcEzpC = fESDZDC->GetZDCP1Energy() / 1000.;
542 Double_t zdcEA = (zdcEznA + zdcEzpA);
543 Double_t zdcEC = (zdcEznC + zdcEzpC);
544 Double_t zdcE = zdcEA + zdcEC;
546 (static_cast<TH1F*>(fHistList->FindObject("fZdcEzdc")))->Fill(zdcE);
547 (static_cast<TH1F*>(fHistList->FindObject("fZdcEzdcAEzdcC")))->Fill(zdcEA, zdcEC);
549 Double_t zdcEzemA = fESDZDC->GetZDCEMEnergy(1) / 1000.;
550 Double_t zdcEzemC = fESDZDC->GetZDCEMEnergy(0) / 1000.;
551 Double_t zdcEzem = zdcEzemA + zdcEzemC;
553 (static_cast<TH2F*>(fHistList->FindObject("fZdcEzemEzdc")))->Fill(zdcEzem, zdcE);
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);
561 // -- ZDC - SPD correlations
563 (static_cast<TH2F*>(fHistList->FindObject("fCorrEzdcSpdOuter")))->Fill(zdcE, fSpdNClustersOuter);
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);
575 //##################################################################################
576 Int_t AliHLTMultiplicityCorrelations::ProcessSPD() {
577 // see header file for class documentation
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);
584 // -- SPD vs TPC correlations
586 (static_cast<TH2F*>(fHistList->FindObject("fCorrSpdOuterTpc")))->Fill(fSpdNClustersOuter,fTpcTracksA);
587 (static_cast<TH2F*>(fHistList->FindObject("fCorrSpdOuterRefTpc")))->Fill(fSpdNClustersOuter,fTpcTracksRef);
593 //##################################################################################
594 Int_t AliHLTMultiplicityCorrelations::ProcessCentrality() {
595 // see header file for class documentation
597 Double_t centV0M = fCentHistV0Mpercentile->GetBinContent(fCentHistV0Mpercentile->FindBin((fVzeroMult)));
598 (static_cast<TH1F*>(fHistList->FindObject("fVzeroCentV0M")))->Fill(centV0M);