4 #include <TClonesArray.h>
11 #include <AliAnalysisTask.h>
12 #include <AliAnalysisManager.h>
14 #include <AliMultiplicity.h>
15 #include <AliESDVertex.h>
16 #include <AliESDEvent.h>
17 #include <AliESDInputHandler.h>
19 #include "AliAnalysisTaskRL.h"
20 #include "AliMCEventHandler.h"
21 #include "AliMCEvent.h"
24 #include "AliTrackReference.h"
26 #include <AliGenEventHeader.h>
27 #include "AliAnalysisTaskSPDdNdEta.h"
30 ClassImp(AliAnalysisTaskSPDdNdEta)
32 //________________________________________________________________________
33 AliAnalysisTaskSPDdNdEta::AliAnalysisTaskSPDdNdEta(const char *name)
34 : AliAnalysisTask(name, "SPDdNdEtaTask"),
42 //Data to be corrected
43 fHistSPDRAWMultvsZ(0),
44 fHistSPDRAWMultvsZTriggEvts(0),
47 //Clusters inner layer and tracklets
48 fHistSPDmultEtacut(0),
51 fHistSPDcl1multEtacutLay1(0),
61 fHistSPDphivsSPDeta(0),
62 fHistSPDcl1phivsSPDcl1eta(0),
68 fHistNcontribSPDvtxvsSPDvtx(0),
69 fHistNcontribSPDvtx3D(0),
70 fHistNcontribSPDvtxZ(0),
71 fHistNcontribSPDvtxall(0),
72 fHistSPDmultvsSPDvtx(0),
75 fHistSPDcl1multvsnFiredChipsLay1(0),
76 fHistSPDmultvsnFiredChipsLay1(0),
77 fHistSPDmultvsnFiredChipsLay2(0),
78 fHistnFiredChipsLay2vsnFiredChipsLay1(0),
79 fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec(0),
81 //Track level correction histograms
87 fHistNonDetectableCorrNum(0),
88 fHistNonDetectableCorrDen(0),
90 fHistTrackTrigVtxCorrNum(0),
91 fHistTrackTrigCorrDen(0),
93 //Event level correction histograms
94 fHistTrigVtxCorrNum(0),
95 fHistTrigVtxCorrDen(0),
100 fHistMCEtavsZTriggMCvtxEvts(0),
101 fHistMCEtavsZTriggESDvtxEvts(0),
106 fHistContributorsvsDeVtx(0),
107 fHistoDetectableNotr(0),
108 fHistoDetectabletr(0),
109 fHistoNonStoppingTracks(0),
110 fHistoDetectedLay1(0),
111 fHistoDetectedLay2(0),
113 fHistoDetectableTRm1(0),
114 fHistoDetectableTR0(0),
115 fHistoDetectableTR1(0),
116 fHistoDetectableTR2(0),
117 fHistoDetectableTR3(0),
118 fHistoDetectableTR4(0),
119 fHistoDetectableTR5(0),
120 fHistoDetectableTR6(0),
127 // Define input and output slots here
129 DefineInput(0, TChain::Class());
130 DefineOutput(0, TList::Class());
133 //________________________________________________________________________
134 AliAnalysisTaskSPDdNdEta::~AliAnalysisTaskSPDdNdEta()
139 // histograms are in the output list and deleted when the output
140 // list is deleted by the TSelector dtor
148 //________________________________________________________________________
149 void AliAnalysisTaskSPDdNdEta::ConnectInputData(Option_t *)
155 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
158 Printf("ERROR: Could not read chain from input slot 0");
160 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
162 Printf("ERROR: Could not get ESDInputHandler");
163 } else fESD = esdH->GetEvent();
165 // disable info messages of AliMCEvent (per event)
166 AliLog::SetClassDebugLevel("AliMCEvent", AliLog::kWarning - AliLog::kDebug + 1);
170 //________________________________________________________________________
171 void AliAnalysisTaskSPDdNdEta::CreateOutputObjects()
178 fHistSPDRAWMultvsZ= new TH2F("fHistSPDRAWMultvsZ", "",200,0.,200.,20,-20.,20.);
179 fHistSPDRAWMultvsZ->Sumw2();
180 fOutput->Add(fHistSPDRAWMultvsZ);
182 fHistSPDRAWMultvsZTriggEvts = new TH2F("fHistSPDRAWMultvsZTriggEvts", "",200,0.,200.,20,-20.,20.);
183 fHistSPDRAWMultvsZTriggEvts->Sumw2();
184 fOutput->Add(fHistSPDRAWMultvsZTriggEvts);
186 fHistSPDRAWEtavsZ = new TH2F("fHistSPDRAWEtavsZ", "Tracklet pseudorapidity distribution", 120, -3.,3.,40,-20.,20.);
187 fHistSPDRAWEtavsZ->GetXaxis()->SetTitle("Pseudorapidity #eta");
188 fHistSPDRAWEtavsZ->GetYaxis()->SetTitle("z_{SPDvtx} [cm]");
189 fHistSPDRAWEtavsZ->Sumw2();
190 fOutput->Add(fHistSPDRAWEtavsZ);
193 fHistSPDmultEtacut = new TH1F("fHistSPDmultEtacut", "Tracklet multiplicity distribution", 200, 0., 200);
194 fHistSPDmultEtacut->GetXaxis()->SetTitle("Reconstructed multiplicity (|#eta|<1.5)");
195 fHistSPDmultEtacut->GetYaxis()->SetTitle("Entries");
196 fOutput->Add(fHistSPDmultEtacut);
198 fHistSPDmult = new TH1F("fHistSPDmult", "Tracklet multiplicity distribution", 200, 0., 200);
199 fHistSPDmult->GetXaxis()->SetTitle("Reconstructed tracklet multiplicity");
200 fHistSPDmult->GetYaxis()->SetTitle("Entries");
201 fOutput->Add(fHistSPDmult);
203 fHistSPDeta = new TH1F("fHistSPDeta", "Tracklet pseudorapidity distribution", 120, -3.,3.);
204 fHistSPDeta->GetXaxis()->SetTitle("Pseudorapidity #eta");
205 fHistSPDeta->GetYaxis()->SetTitle("Entries");
206 fHistSPDeta->SetLineColor(kGreen);
207 fHistSPDeta->SetLineWidth(3);
208 fHistSPDeta->Sumw2();
209 fOutput->Add(fHistSPDeta);
211 fHistSPDcl1multEtacutLay1 = new TH1F("fHistSPDcl1multEtacutLay1", "Cluster multiplicity (inner layer)", 200, 0., 200);
212 fHistSPDcl1multEtacutLay1->GetXaxis()->SetTitle("Cluster multiplicity lay1 (|#eta|<2.)");
213 fHistSPDcl1multEtacutLay1->GetYaxis()->SetTitle("Entries");
214 fOutput->Add(fHistSPDcl1multEtacutLay1);
216 fHistSPDcl1mult = new TH1F("fHistSPDcl1mult", "Cluster multiplicity (inner layer)", 200, 0., 200);
217 fHistSPDcl1mult->GetXaxis()->SetTitle("Cluster multiplicity lay1");
218 fHistSPDcl1mult->GetYaxis()->SetTitle("Entries");
219 fOutput->Add(fHistSPDcl1mult);
221 fHistSPDcl1eta = new TH1F("fHistSPDcl1eta", "Cluster pseudorapidity (inner layer)", 120, -3.,3.);
222 fHistSPDcl1eta->GetXaxis()->SetTitle("Pseudorapidity #eta");
223 fHistSPDcl1eta->GetYaxis()->SetTitle("Entries");
224 fHistSPDcl1eta->Sumw2();
225 fOutput->Add(fHistSPDcl1eta);
227 fHistSPDphi = new TH1F("fHistSPDphi", "Tracklet #phi distribution", 360, 0.,2*TMath::Pi());
228 fHistSPDphi->GetXaxis()->SetTitle("#varphi [rad]");
229 fHistSPDphi->GetYaxis()->SetTitle("Entries");
230 fOutput->Add(fHistSPDphi);
232 fHistSPDcl1phi= new TH1F("fHistSPDcl1phi", "Cluster #phi (inner layer) ", 360, 0.,2*TMath::Pi());
233 fHistSPDcl1phi->GetXaxis()->SetTitle("#varphi [rad]");
234 fHistSPDcl1phi->GetYaxis()->SetTitle("Entries");
235 fOutput->Add(fHistSPDcl1phi);
237 fHistSPDtheta = new TH1F("fHistSPDtheta", "Tracklet #theta distribution", 360, 0.,2*TMath::Pi());
238 fHistSPDtheta->GetXaxis()->SetTitle("#theta [rad]");
239 fHistSPDtheta->GetYaxis()->SetTitle("Entries");
240 fOutput->Add(fHistSPDtheta);
242 fHistSPDcl1theta = new TH1F("fHistSPDcl1theta", "Cluster #theta (inner layer)", 360, 0.,2*TMath::Pi());
243 fHistSPDcl1theta->GetXaxis()->SetTitle("#theta [rad]");
244 fHistSPDcl1theta->GetYaxis()->SetTitle("Entries");
245 fOutput->Add(fHistSPDcl1theta);
247 fHistSPDdePhi= new TH1F("fHistSPDdePhi", "Tracklet #Delta#varphi distribution",400,-0.1,.1);
248 fHistSPDdePhi->GetXaxis()->SetTitle("#Delta#varphi [rad]");
249 fHistSPDdePhi->GetYaxis()->SetTitle("z_{MCvtx}");
250 fOutput->Add(fHistSPDdePhi);
252 fHistSPDdePhiZ= new TH1F("fHistSPDdePhiZ", "Tracklet #Delta#varphi distribution",400,-0.1,.1);
253 fOutput->Add(fHistSPDdePhiZ);
255 fHistSPDdePhi3D= new TH1F("fHistSPDdePhi3D", "Tracklet #Delta#varphi distribution",400,-0.1,.1);
256 fOutput->Add(fHistSPDdePhi3D);
258 fHistSPDphivsSPDeta= new TH2F("fHistSPDphivsSPDeta", "Tracklets - #varphi vs #eta",120,-3.,3,360,0.,2*TMath::Pi());
259 fHistSPDphivsSPDeta->GetXaxis()->SetTitle("Pseudorapidity #eta");
260 fHistSPDphivsSPDeta->GetYaxis()->SetTitle("#varphi [rad]");
261 fOutput->Add(fHistSPDphivsSPDeta);
263 fHistSPDcl1phivsSPDcl1eta= new TH2F("fHistSPDcl1phivsSPDcl1eta", "Clusters layer1 - #varphi vs #eta",120,-3.,3,360,0.,2*TMath::Pi());
264 fHistSPDcl1phivsSPDcl1eta->GetXaxis()->SetTitle("Pseudorapidity #eta");
265 fHistSPDcl1phivsSPDcl1eta->GetYaxis()->SetTitle("#varphi [rad]");
266 fOutput->Add(fHistSPDcl1phivsSPDcl1eta);
269 fHistSPDvtx = new TH1F("fHistSPDvtx", "SPD vertex distribution - all events",20,-20.,20.);
270 fHistSPDvtx->GetXaxis()->SetTitle("z_{SPDvtx} [cm]");
271 fHistSPDvtx->GetYaxis()->SetTitle("Entries");
272 fOutput->Add(fHistSPDvtx);
274 fHistSPDvtx3D = new TH3F("fHistSPDvtx3D", "SPD vertex distribution",100,-5.,5.,100,-5.,5.,400,-20.,20.);
275 fOutput->Add(fHistSPDvtx3D);
277 fHistSPDvtxZ = new TH3F("fHistSPDvtxZ", "SPD vertex distribution",100,-5.,5.,100,-5.,5.,400,-20.,20.);
278 fOutput->Add(fHistSPDvtxZ);
280 fHistNcontribSPDvtxvsSPDvtx= new TH2F("fHistNcontribSPDvtxvsSPDvtx", " ",100,-50.,50.,202,-2.,200.);
281 fHistNcontribSPDvtxvsSPDvtx->GetXaxis()->SetTitle("z_{SPDvtx} [cm]");
282 fHistNcontribSPDvtxvsSPDvtx->GetYaxis()->SetTitle("# contributors");
283 fOutput->Add(fHistNcontribSPDvtxvsSPDvtx);
285 fHistNcontribSPDvtx3D= new TH1F("fHistNcontribSPDvtx_3D", "SPD vtx 3D",202,-2.,200.);
286 fHistNcontribSPDvtx3D->GetXaxis()->SetTitle("# contributors");
287 fHistNcontribSPDvtx3D->GetYaxis()->SetTitle("Entries");
288 fOutput->Add(fHistNcontribSPDvtx3D);
290 fHistNcontribSPDvtxZ= new TH1F("fHistNcontribSPDvtxZ", "SPD vtx Z",202,-2.,200.);
291 fHistNcontribSPDvtxZ->GetXaxis()->SetTitle("# contributors");
292 fHistNcontribSPDvtxZ->GetYaxis()->SetTitle("Entries");
293 fOutput->Add(fHistNcontribSPDvtxZ);
295 fHistNcontribSPDvtxall= new TH1F("fHistNcontribSPDvtxall", "SPD vtx - all events",202,-2.,200.);
296 fHistNcontribSPDvtxall->GetXaxis()->SetTitle("# contributors");
297 fHistNcontribSPDvtxall->GetYaxis()->SetTitle("Entries");
298 fOutput->Add(fHistNcontribSPDvtxall);
300 fHistSPDmultvsSPDvtx= new TH2F("fHistSPDmultvsSPDvtx", "",20,-20.,20.,200,0.,200.);
301 fHistSPDmultvsSPDvtx->GetXaxis()->SetTitle("z_{recvtx} [cm]");
302 fHistSPDmultvsSPDvtx->GetYaxis()->SetTitle("Reconstructed multiplicity");
303 fOutput->Add(fHistSPDmultvsSPDvtx);
305 fHistSPDcl1multvsnFiredChipsLay1 = new TH2F("fHistSPDcl1multvsnFiredChipsLay1", "",200,0.,200.,200,0.,200.);
306 fHistSPDcl1multvsnFiredChipsLay1->GetXaxis()->SetTitle("# fired chips lay1");
307 fHistSPDcl1multvsnFiredChipsLay1->GetYaxis()->SetTitle("Cluster lay1 multiplicity");
308 fOutput->Add(fHistSPDcl1multvsnFiredChipsLay1);
310 fHistSPDmultvsnFiredChipsLay1 = new TH2F("fHistSPDmultvsnFiredChipsLay1","",200,0.,200.,200,0.,200.);
311 fHistSPDmultvsnFiredChipsLay1->GetXaxis()->SetTitle("# fired chips lay1");
312 fHistSPDmultvsnFiredChipsLay1->GetYaxis()->SetTitle("Tracklet multiplicity");
313 fOutput->Add(fHistSPDmultvsnFiredChipsLay1);
315 fHistSPDmultvsnFiredChipsLay2 = new TH2F("fHistSPDmultvsnFiredChipsLay2","",200,0.,200.,200,0.,200.);
316 fHistSPDmultvsnFiredChipsLay2->GetXaxis()->SetTitle("# fired chips lay2");
317 fHistSPDmultvsnFiredChipsLay2->GetYaxis()->SetTitle("Tracklet multiplicity");
318 fOutput->Add(fHistSPDmultvsnFiredChipsLay2);
320 fHistnFiredChipsLay2vsnFiredChipsLay1 = new TH2F("fHistnFiredChipsLay2vsnFiredChipsLay1","",200,0.,200.,200,0.,200.);
321 fHistnFiredChipsLay2vsnFiredChipsLay1->GetXaxis()->SetTitle("# fired chips lay1");
322 fHistnFiredChipsLay2vsnFiredChipsLay1->GetYaxis()->SetTitle("# fired chip lay2");
323 fOutput->Add(fHistnFiredChipsLay2vsnFiredChipsLay1);
325 fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec = new TH2F("fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec","",200,0.,200.,200,0.,200.);
326 fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec->GetXaxis()->SetTitle("# fired chips lay1");
327 fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec->GetYaxis()->SetTitle("# fired chip lay2");
328 fOutput->Add(fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec);
332 fHistBkgCorrNum = new TH2F("fHistBkgCorrNum","",60,-3.00,3.00,20,-20.,20.);
333 fHistBkgCorrNum->Sumw2();
334 fOutput->Add(fHistBkgCorrNum);
336 fHistBkgCorrDen = new TH2F("fHistBkgCorrDen","",60,-3.00,3.00,20,-20.,20.);
337 fHistBkgCorrDen->Sumw2();
338 fOutput->Add(fHistBkgCorrDen);
340 fHistAlgEffNum = new TH2F("fHistAlgEffNum","",60,-3.00,3.00,20,-20.,20.);
341 fHistAlgEffNum->Sumw2();
342 fOutput->Add(fHistAlgEffNum);
344 fHistNonDetectableCorrNum = new TH2F("fHistNonDetectableCorrNum","",60,-3.00,3.00,20,-20.,20.);
345 fHistNonDetectableCorrNum->Sumw2();
346 fOutput->Add(fHistNonDetectableCorrNum);
348 fHistNonDetectableCorrDen = new TH2F("fHistNonDetectableCorrDen","",60,-3.00,3.00,20,-20.,20.);
349 fHistNonDetectableCorrDen->Sumw2();
350 fOutput->Add(fHistNonDetectableCorrDen);
352 fHistTrackTrigVtxCorrNum = new TH2F("fHistTrackTrigVtxCorrNum","",60,-3.00,3.00,20,-20.,20.);
353 fHistTrackTrigVtxCorrNum->Sumw2();
354 fOutput->Add(fHistTrackTrigVtxCorrNum);
356 fHistTrackTrigCorrDen = new TH2F("fHistTrackTrigCorrDen","",60,-3.00,3.00,20,-20.,20.);
357 fHistTrackTrigCorrDen->Sumw2();
358 fOutput->Add(fHistTrackTrigCorrDen);
360 //Event level correction histograms
361 fHistTrigVtxCorrNum = new TH2F("fHistTrigVtxCorrNum","",200,0.,200.,20,-20.,20.);
362 fHistTrigVtxCorrNum->Sumw2();
363 fOutput->Add(fHistTrigVtxCorrNum);
365 fHistTrigVtxCorrDen = new TH2F("fHistTrigVtxCorrDen","",200,0.,200.,20,-20.,20.);
366 fHistTrigVtxCorrDen->Sumw2();
367 fOutput->Add(fHistTrigVtxCorrDen);
369 fHistTrigCorrDen = new TH2F("fHistTrigCorrDen","",200,0.,200.,20,-20.,20.);
370 fHistTrigCorrDen->Sumw2();
371 fOutput->Add(fHistTrigCorrDen);
374 fHistMCEtavsZTriggMCvtxEvts = new TH2F("fHistMCEtavsZTriggMCvtxEvts","Generated pseudorapidity distribution",60,-3.00,3.00,20,-20.,20.);
375 fHistMCEtavsZTriggMCvtxEvts->GetXaxis()->SetTitle("Pseudorapidity #eta");
376 fHistMCEtavsZTriggMCvtxEvts->GetYaxis()->SetTitle("Entries");
377 fHistMCEtavsZTriggMCvtxEvts->Sumw2();
378 fOutput->Add(fHistMCEtavsZTriggMCvtxEvts);
380 fHistMCEtavsZTriggESDvtxEvts = new TH2F("fHistMCEtavsZTriggESDvtxEvts","Generated pseudorapidity distribution",60,-3.00,3.00,20,-20.,20.);
381 fHistMCEtavsZTriggESDvtxEvts->GetXaxis()->SetTitle("Pseudorapidity #eta");
382 fHistMCEtavsZTriggESDvtxEvts->GetYaxis()->SetTitle("Entries");
383 fHistMCEtavsZTriggESDvtxEvts->Sumw2();
384 fOutput->Add(fHistMCEtavsZTriggESDvtxEvts);
386 fHistMCEtavsZ = new TH2F("fHistMCEtavsZ","Generated pseudorapidity distribution",60,-3.00,3.00,20,-20.,20.);
387 fHistMCEtavsZ->GetXaxis()->SetTitle("Pseudorapidity #eta");
388 fHistMCEtavsZ->GetYaxis()->SetTitle("Entries");
389 fHistMCEtavsZ->Sumw2();
390 fOutput->Add(fHistMCEtavsZ);
392 fHistTRradius = new TH1F("fHistTRradius","ITS track reference rad",200,0.,10.);
393 fOutput->Add(fHistTRradius);
395 fHistContributorsvsDeVtx = new TH2F("fHistContributorsvsDeVtx","",200,-20.,20.,22,-2.,20.);
396 fOutput->Add(fHistContributorsvsDeVtx);
398 fHistoDetectableNotr = new TH3F("fHistoDetectableNotr","",60,-3.00,3.00,20,-20.,20.,100,0.,10.);
399 fOutput->Add(fHistoDetectableNotr);
401 fHistoDetectabletr = new TH2F("fHistoDetectabletr","",60,-3.00,3.00,20,-20.,20.);
402 fOutput->Add(fHistoDetectabletr);
404 fHistoNonStoppingTracks = new TH2F("fHistoNonStoppingTracks","",60,-3.00,3.00,20,-20.,20.);
405 fOutput->Add(fHistoNonStoppingTracks);
407 fHistoDetectedLay1 = new TH2F("fHistoDetectedLay1","",60,-3.00,3.00,20,-20.,20.);
408 fOutput->Add(fHistoDetectedLay1);
410 fHistoDetectedLay2 = new TH2F("fHistoDetectedLay2","",60,-3.00,3.00,20,-20.,20.);
411 fOutput->Add(fHistoDetectedLay2);
413 fHistoPt = new TH1F("fHistoPt","",100,.0,10.);
414 fOutput->Add(fHistoPt);
416 fHistoDetectableTRm1 = new TH2F("fHistoDetectableTRm1","",60,-3.00,3.00,20,-20.,20.);
417 fOutput->Add(fHistoDetectableTRm1);
419 fHistoDetectableTR0 = new TH2F("fHistoDetectableTR0","",60,-3.00,3.00,20,-20.,20.);
420 fOutput->Add(fHistoDetectableTR0);
422 fHistoDetectableTR1 = new TH2F("fHistoDetectableTR1","",60,-3.00,3.00,20,-20.,20.);
423 fOutput->Add(fHistoDetectableTR1);
425 fHistoDetectableTR2 = new TH2F("fHistoDetectableTR2","",60,-3.00,3.00,20,-20.,20.);
426 fOutput->Add(fHistoDetectableTR2);
428 fHistoDetectableTR3 = new TH2F("fHistoDetectableTR3","",60,-3.00,3.00,20,-20.,20.);
429 fOutput->Add(fHistoDetectableTR3);
431 fHistoDetectableTR4 = new TH2F("fHistoDetectableTR4","",60,-3.00,3.00,20,-20.,20.);
432 fOutput->Add(fHistoDetectableTR4);
434 fHistoDetectableTR5 = new TH2F("fHistoDetectableTR5","",60,-3.00,3.00,20,-20.,20.);
435 fOutput->Add(fHistoDetectableTR5);
437 fHistoDetectableTR6 = new TH2F("fHistoDetectableTR6","",60,-3.00,3.00,20,-20.,20.);
438 fOutput->Add(fHistoDetectableTR6);
440 fHistoRTRm1 = new TH1F("fHistoRTRm1","",10000,0.,5000);
441 fOutput->Add(fHistoRTRm1);
446 //________________________________________________________________________
447 void AliAnalysisTaskSPDdNdEta::Exec(Option_t *)
450 // Called for each event
454 Printf("ERROR: fESD not available");
459 const AliESDVertex* vtxESD = fESD->GetVertex();
461 vtxESD->GetXYZ(esdvtx);
462 Int_t nContrib = vtxESD->GetNContributors();
464 const AliMultiplicity* multESD = fESD->GetMultiplicity();
466 // Loading tracklets...
468 multSPD = multESD->GetNumberOfTracklets();
471 ULong64_t triggerMask;
472 ULong64_t spdFO = (1 << 14);
473 ULong64_t v0left = (1 << 11);
474 ULong64_t v0right = (1 << 12);
476 triggerMask=fESD->GetTriggerMask();
478 Bool_t eventTriggered = kFALSE;
480 if (fTrigger==0) eventTriggered = kTRUE;
482 if (fTrigger==1) eventTriggered = triggerMask&spdFO || ((triggerMask&v0left) || (triggerMask&v0right));
484 if (fTrigger==2) eventTriggered = triggerMask&spdFO && ((triggerMask&v0left) || (triggerMask&v0right));
486 PostData(0, fOutput);
488 AliMultiplicity * mult1 = (AliMultiplicity*)multESD;
489 Short_t nFiredChipsLay1 = mult1->GetNumberOfFiredChips(0);
490 Short_t nFiredChipsLay2 = mult1->GetNumberOfFiredChips(1);
491 Int_t multSPDEtacut = 0;
492 Int_t multSPDcl1 = 0;
493 Int_t nSingleCl1 = 0;
494 Int_t multSPDcl1EtacutLay1 = 0;
495 nSingleCl1 = multESD->GetNumberOfSingleClusters();
496 multSPDcl1 = nSingleCl1 + multSPD;
497 Float_t* recEtaSPDcl1 = new Float_t[multSPD+nSingleCl1];
499 // Printf("There are %d tracklets in this event", multSPD);
500 if (esdvtx[2]!=0.&&eventTriggered&&multSPD!=0) {
502 fHistSPDvtx->Fill(esdvtx[2]);
503 if (strcmp(vtxESD->GetTitle(),"vertexer: 3D") == 0) fHistSPDvtx3D->Fill(esdvtx[0],esdvtx[1],esdvtx[2]);
504 else fHistSPDvtxZ->Fill(esdvtx[0],esdvtx[1],esdvtx[2]);
506 for (Int_t itracklet=0; itracklet<multSPD; ++itracklet) {
507 Float_t thetaTr= multESD->GetTheta(itracklet);
508 Float_t phiTr= multESD->GetPhi(itracklet);
509 Float_t dePhiTr= multESD->GetDeltaPhi(itracklet);
510 Float_t recEtaSPD =multESD->GetEta(itracklet);
511 recEtaSPDcl1[itracklet] = recEtaSPD;
513 fHistSPDeta->Fill(recEtaSPD);
514 fHistSPDRAWEtavsZ->Fill(recEtaSPD,esdvtx[2]);
515 fHistSPDcl1eta->Fill(recEtaSPD);
516 fHistSPDphi->Fill(phiTr);
517 fHistSPDcl1phi->Fill(phiTr);
518 fHistSPDtheta->Fill(thetaTr);
519 fHistSPDcl1theta->Fill(thetaTr);
520 fHistSPDdePhi->Fill(dePhiTr);
521 if (strcmp(vtxESD->GetTitle(),"vertexer: Z") == 0) fHistSPDdePhiZ->Fill(dePhiTr);
522 if (strcmp(vtxESD->GetTitle(),"vertexer: 3D") == 0) fHistSPDdePhi3D->Fill(dePhiTr);
523 fHistSPDphivsSPDeta->Fill(recEtaSPD,phiTr);
524 fHistSPDcl1phivsSPDcl1eta->Fill(recEtaSPD,phiTr);
526 // calculate multiplicity in etacut
527 if (TMath::Abs(recEtaSPD)<1.5) multSPDEtacut++;
528 if (TMath::Abs(recEtaSPD)<2.) multSPDcl1EtacutLay1++;
531 for (Int_t iCl1=0; iCl1<nSingleCl1; ++iCl1) {
532 Float_t thetaSingleCl1 = multESD->GetThetaSingle(iCl1);
534 Float_t etaSingleCl1 = -TMath::Log(TMath::Tan(thetaSingleCl1/2.));
535 Float_t phiSingleCl1 = multESD->GetPhiSingle(iCl1);
536 recEtaSPDcl1[iCl1+multSPD] = etaSingleCl1;
538 fHistSPDcl1eta->Fill(etaSingleCl1);
539 fHistSPDcl1phi->Fill(phiSingleCl1);
540 fHistSPDcl1theta->Fill(thetaSingleCl1);
541 fHistSPDcl1phivsSPDcl1eta->Fill(etaSingleCl1,phiSingleCl1);
542 if (TMath::Abs(etaSingleCl1)<2.) multSPDcl1EtacutLay1++;
545 fHistSPDmultEtacut->Fill(multSPDEtacut);
546 fHistSPDmult->Fill(multSPD);
547 fHistSPDcl1multEtacutLay1->Fill(multSPDcl1EtacutLay1);
548 if (strcmp(vtxESD->GetTitle(),"vertexer: 3D") == 0) {
549 fHistNcontribSPDvtx3D->Fill(nContrib);
551 if (strcmp(vtxESD->GetTitle(),"vertexer: Z") == 0) {
552 fHistNcontribSPDvtxZ->Fill(nContrib);
554 fHistSPDmultvsSPDvtx->Fill(esdvtx[2],multSPD);
555 fHistNcontribSPDvtxvsSPDvtx->Fill(esdvtx[2],nContrib);
556 fHistSPDRAWMultvsZ->Fill(multSPD,esdvtx[2]);
557 fHistSPDmultvsnFiredChipsLay1->Fill(nFiredChipsLay1,multSPD);
558 fHistSPDmultvsnFiredChipsLay2->Fill(nFiredChipsLay2,multSPD);
562 if (eventTriggered) fHistSPDRAWMultvsZTriggEvts->Fill(multSPD,esdvtx[2]);
564 fHistSPDcl1mult->Fill(multSPDcl1);
565 fHistNcontribSPDvtxall->Fill(nContrib);
567 fHistSPDcl1multvsnFiredChipsLay1->Fill(nFiredChipsLay1,multSPDcl1);
569 fHistnFiredChipsLay2vsnFiredChipsLay1->Fill(nFiredChipsLay1,nFiredChipsLay2);
572 fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec->Fill(nFiredChipsLay1,nFiredChipsLay2);
575 delete[] recEtaSPDcl1;
578 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
580 Printf("ERROR: Could not retrieve MC event handler");
584 AliMCEvent* mcEvent = eventHandler->MCEvent();
586 Printf("ERROR: Could not retrieve MC event");
590 AliStack* stack = mcEvent->Stack();
592 AliDebug(AliLog::kError, "Stack not available");
596 AliHeader* header = mcEvent->Header();
598 AliDebug(AliLog::kError, "Header not available");
601 AliGenEventHeader* genHeader = header->GenEventHeader();
604 genHeader->PrimaryVertex(vtxMC);
608 Int_t multMCCharged = 0;
609 Int_t nMCPart = stack->GetNprimary();
610 Float_t* eta_gen = new Float_t[nMCPart];
611 Float_t* pt_gen = new Float_t[nMCPart];
612 Int_t* stackIndexOfPrimaryParts = new Int_t[nMCPart];
613 Int_t* reconstructedPrimaryPart = new Int_t[nMCPart];
614 Int_t* detectedPrimaryPartLay1 = new Int_t[nMCPart];
615 Int_t* detectedPrimaryPartLay2 = new Int_t[nMCPart];
616 Int_t* detectablePrimaryPart = new Int_t[nMCPart];
617 Bool_t* stoppingPartAtLay2 = new Bool_t[nMCPart];
619 // Loading track references...
620 Float_t rminL1 = 3.4;
621 Float_t rmaxL1 = 4.4;
622 Float_t rminL2 = 6.9;
623 Float_t rmaxL2 = 7.9;
625 Int_t* nTrackRefPerPrim = new Int_t[nMCPart];
626 Float_t trackRefPerPrim[nMCPart][500][2];
628 TTree* tRefTree = eventHandler->TreeTR();
630 AliTrackReference *tref=0x0;
631 mcEvent->ConnectTreeTR(tRefTree);
633 // Loop over MC particles
634 for (Int_t imc=0; imc<nMCPart; imc++) {
635 TParticle* part = stack->Particle(imc);
636 Bool_t IsPrimary = stack->IsPhysicalPrimary(imc);
637 if (!IsPrimary) continue;
638 TParticlePDG* pdgPart = part->GetPDG();
639 if (TMath::Abs(pdgPart->Charge())!=3) continue;
640 Float_t theta = part->Theta();
641 if (theta==0 || theta==TMath::Pi()) continue;
642 Float_t eta = part->Eta();
643 Float_t pt = part->Pt();
644 eta_gen[multMCCharged] = eta;
645 pt_gen[multMCCharged] = pt;
646 stackIndexOfPrimaryParts[multMCCharged] = imc;
648 reconstructedPrimaryPart[multMCCharged]=kFALSE;
649 detectedPrimaryPartLay1[multMCCharged]=kFALSE;
650 detectedPrimaryPartLay2[multMCCharged]=kFALSE;
651 detectablePrimaryPart[multMCCharged]=kFALSE;
652 stoppingPartAtLay2[multMCCharged]=kFALSE;
654 fHistoPt->Fill(pt_gen[multMCCharged]);
656 AliMCParticle* mcpart = mcEvent->GetTrack(imc);
657 Int_t nref = mcpart->GetNumberOfTrackReferences();
658 nTrackRefPerPrim[multMCCharged]=nref;
660 Bool_t nonStoppingPart=kFALSE;
662 // Loop over all the track refs of the track
663 for (Int_t iref=0;iref<nref;iref++) {
664 tref = mcpart->GetTrackReference(iref);
666 trackRefPerPrim[multMCCharged][iref][0]=tref->R();
667 trackRefPerPrim[multMCCharged][iref][1]=(Float_t)tref->DetectorId();
668 if (trackRefPerPrim[multMCCharged][iref][1]==0) fHistTRradius->Fill(trackRefPerPrim[multMCCharged][iref][0]); // ITS only
669 // Printf("det ID: %f R: %f", trackRefPerPrim[multMCCharged][iref][1], trackRefPerPrim[multMCCharged][iref][0]);
670 if (iref==nref-1&&trackRefPerPrim[multMCCharged][iref][1]!=-1) {
671 nonStoppingPart=kTRUE;
675 if (nref==0||nonStoppingPart) {
676 detectablePrimaryPart[multMCCharged]=kTRUE; //0. if prim doesn't disappear it is detectable
677 if (nonStoppingPart) fHistoNonStoppingTracks->Fill(eta_gen[multMCCharged],vtxMC[2]);
678 if (nTrackRefPerPrim[multMCCharged]==0) fHistoDetectableNotr->Fill(eta_gen[multMCCharged],vtxMC[2],pt_gen[multMCCharged]);
682 } //end of MC particle loop
685 if (multSPD==0 && eventTriggered && esdvtx[2]!=0.) {
686 fHistContributorsvsDeVtx->Fill(esdvtx[2]-vtxMC[2],nContrib);
690 if (eventTriggered && esdvtx[2]!=0. && multSPD!=0) {
692 for (Int_t imc=0; imc<multMCCharged; imc++) {
694 // Find out detected prims on each layer
695 for (Int_t iref=0; iref<nTrackRefPerPrim[imc]; iref++) {
696 if (trackRefPerPrim[imc][iref][0]>rminL1&&trackRefPerPrim[imc][iref][0]<rmaxL1&&trackRefPerPrim[imc][iref][1]==0) {
697 detectedPrimaryPartLay1[imc] = kTRUE;
698 fHistoDetectedLay1->Fill(eta_gen[imc],vtxMC[2]);
702 for (Int_t iref=0; iref<nTrackRefPerPrim[imc]; iref++) {
703 if (trackRefPerPrim[imc][iref][0]>rminL2&&trackRefPerPrim[imc][iref][0]<rmaxL2) {
704 if (trackRefPerPrim[imc][iref][1]==0) {
705 detectedPrimaryPartLay2[imc] = kTRUE;
706 fHistoDetectedLay2->Fill(eta_gen[imc],vtxMC[2]);
707 detectablePrimaryPart[imc] = kTRUE; //1. if prim is detected it is also detectable
710 if (trackRefPerPrim[imc][iref][1]==-1) {
711 stoppingPartAtLay2[imc]=kTRUE;
716 // Find out SPD detectable prims
717 if (!detectablePrimaryPart[imc]) {
718 for (Int_t iref=0; iref<nTrackRefPerPrim[imc]; iref++) {
719 if (trackRefPerPrim[imc][iref][0]>rmaxL2) {
720 detectablePrimaryPart[imc] = kTRUE; //2. if it has crossed layer 2 it is detectable
721 fHistoDetectabletr->Fill(eta_gen[imc],vtxMC[2]);
722 if (trackRefPerPrim[imc][iref][1]==-1) {
723 fHistoDetectableTRm1->Fill(eta_gen[imc],vtxMC[2]);
724 fHistoRTRm1->Fill(trackRefPerPrim[imc][iref][0]);
726 else if (trackRefPerPrim[imc][iref][1]==0) fHistoDetectableTR0->Fill(eta_gen[imc],vtxMC[2]);
727 else if (trackRefPerPrim[imc][iref][1]==1) fHistoDetectableTR1->Fill(eta_gen[imc],vtxMC[2]);
728 else if (trackRefPerPrim[imc][iref][1]==2) fHistoDetectableTR2->Fill(eta_gen[imc],vtxMC[2]);
729 else if (trackRefPerPrim[imc][iref][1]==3) fHistoDetectableTR3->Fill(eta_gen[imc],vtxMC[2]);
730 else if (trackRefPerPrim[imc][iref][1]==4) fHistoDetectableTR4->Fill(eta_gen[imc],vtxMC[2]);
731 else if (trackRefPerPrim[imc][iref][1]==5) fHistoDetectableTR5->Fill(eta_gen[imc],vtxMC[2]);
732 else if (trackRefPerPrim[imc][iref][1]==6) fHistoDetectableTR6->Fill(eta_gen[imc],vtxMC[2]);
737 // Printf("Primary nr %d", stackIndexOfPrimaryParts[imc]);
740 for (Int_t itracklet=0; itracklet<multSPD; ++itracklet) {
742 Int_t labL1 = multESD->GetLabel(itracklet,0);
743 Int_t labL2 = multESD->GetLabel(itracklet,1);
745 fHistBkgCorrDen->Fill(multESD->GetEta(itracklet),esdvtx[2]);
748 for (Int_t imc=0; imc<multMCCharged; imc++) {
749 if (labL1==stackIndexOfPrimaryParts[imc]) {
750 if (detectedPrimaryPartLay1[imc]&&detectedPrimaryPartLay2[imc]) {
751 reconstructedPrimaryPart[imc]=kTRUE;
759 for (Int_t imc=0; imc<multMCCharged; imc++) {
760 if (reconstructedPrimaryPart[imc]) {
761 fHistBkgCorrNum->Fill(eta_gen[imc],vtxMC[2]);
763 if (detectedPrimaryPartLay1[imc]&&detectedPrimaryPartLay2[imc]) fHistAlgEffNum->Fill(eta_gen[imc],vtxMC[2]);
765 fHistNonDetectableCorrNum->Fill(eta_gen[imc],vtxMC[2]);
767 if (detectablePrimaryPart[imc]) fHistNonDetectableCorrDen->Fill(eta_gen[imc],vtxMC[2]);
769 fHistMCEtavsZTriggESDvtxEvts->Fill(eta_gen[imc],esdvtx[2]);
770 fHistMCEtavsZTriggMCvtxEvts->Fill(eta_gen[imc],vtxMC[2]);
773 fHistTrigVtxCorrDen->Fill(multSPD,vtxMC[2]);
777 if (eventTriggered) {
778 fHistTrigCorrDen->Fill(multSPD,vtxMC[2]);
779 for (Int_t imc=0; imc<multMCCharged; imc++) {
780 fHistTrackTrigCorrDen->Fill(eta_gen[imc],vtxMC[2]);
784 for (Int_t imc=0; imc<multMCCharged; imc++) {
785 fHistTrackTrigVtxCorrNum->Fill(eta_gen[imc],vtxMC[2]);
786 fHistMCEtavsZ->Fill(eta_gen[imc],vtxMC[2]);
789 fHistTrigVtxCorrNum->Fill(multSPD,vtxMC[2]);
793 delete[] nTrackRefPerPrim;
794 delete[] stackIndexOfPrimaryParts;
795 delete[] reconstructedPrimaryPart;
796 delete[] detectedPrimaryPartLay1;
797 delete[] detectedPrimaryPartLay2;
798 delete[] detectablePrimaryPart;
799 delete[] stoppingPartAtLay2;
803 //________________________________________________________________________
804 void AliAnalysisTaskSPDdNdEta::Terminate(Option_t *)
806 // Called once at the end of the query
807 fOutput = dynamic_cast<TList*> (GetOutputData(0));
808 if (!fOutput) { Printf("ERROR: fOutput not available");
812 fHistSPDRAWMultvsZ= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDRAWMultvsZ"));
813 fHistSPDRAWMultvsZTriggEvts = dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDRAWMultvsZTriggEvts"));
814 fHistSPDRAWEtavsZ = dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDRAWEtavsZ"));
816 fHistSPDmultEtacut = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDmultEtacut"));
817 fHistSPDmult = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDmult"));
818 fHistSPDeta = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDeta"));
819 fHistSPDcl1multEtacutLay1 = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDcl1multEtacutLay1"));
820 fHistSPDcl1mult = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDcl1mult"));
821 fHistSPDcl1eta = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDcl1eta"));
822 fHistSPDphi = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDphi"));
823 fHistSPDcl1phi = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDcl1phi"));
824 fHistSPDtheta = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDtheta"));
825 fHistSPDcl1theta = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDcl1theta"));
826 fHistSPDdePhi = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDdePhi"));
827 fHistSPDdePhiZ = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDdePhiZ"));
828 fHistSPDdePhi3D = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDdePhi3D"));
829 fHistSPDphivsSPDeta= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDphivsSPDeta"));
830 fHistSPDcl1phivsSPDcl1eta= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDcl1phivsSPDcl1eta"));
832 fHistSPDvtx = dynamic_cast<TH1F*> (fOutput->FindObject("fHistSPDvtx"));
833 fHistSPDvtx3D = dynamic_cast<TH3F*> (fOutput->FindObject("fHistSPDvtx3D"));
834 fHistSPDvtxZ = dynamic_cast<TH3F*> (fOutput->FindObject("fHistSPDvtxZ"));
835 fHistNcontribSPDvtxvsSPDvtx = dynamic_cast<TH2F*> (fOutput->FindObject("fHistNcontribSPDvtxvsSPDvtx"));
836 fHistNcontribSPDvtx3D = dynamic_cast<TH1F*> (fOutput->FindObject("fHistNcontribSPDvtx3D"));
837 fHistNcontribSPDvtxZ = dynamic_cast<TH1F*> (fOutput->FindObject("fHistNcontribSPDvtxZ"));
838 fHistNcontribSPDvtxall = dynamic_cast<TH1F*> (fOutput->FindObject("fHistNcontribSPDvtxall"));
839 fHistSPDmultvsSPDvtx= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDmultvsSPDvtx"));
841 fHistSPDcl1multvsnFiredChipsLay1= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDcl1multvsnFiredChipsLay1"));
842 fHistSPDmultvsnFiredChipsLay1= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDmultvsnFiredChipsLay1"));
843 fHistSPDmultvsnFiredChipsLay2= dynamic_cast<TH2F*> (fOutput->FindObject("fHistSPDmultvsnFiredChipsLay2"));
844 fHistnFiredChipsLay2vsnFiredChipsLay1= dynamic_cast<TH2F*> (fOutput->FindObject("fHistnFiredChipsLay2vsnFiredChipsLay1"));
845 fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec= dynamic_cast<TH2F*> (fOutput->FindObject("fHistnFiredChipsLay2vsnFiredChipsLay1novtxrec"));
848 fHistBkgCorrNum = dynamic_cast<TH2F*> (fOutput->FindObject("fHistBkgCorrNum"));
849 fHistBkgCorrDen = dynamic_cast<TH2F*> (fOutput->FindObject("fHistBkgCorrDen"));
851 fHistAlgEffNum = dynamic_cast<TH2F*> (fOutput->FindObject("fHistAlgEffNum"));
853 fHistNonDetectableCorrNum = dynamic_cast<TH2F*> (fOutput->FindObject("fHistNonDetectableCorrNum"));
854 fHistNonDetectableCorrDen = dynamic_cast<TH2F*> (fOutput->FindObject("fHistNonDetectableCorrDen"));
856 fHistTrackTrigVtxCorrNum = dynamic_cast<TH2F*> (fOutput->FindObject("fHistTrackTrigVtxCorrNum"));
858 fHistTrackTrigCorrDen = dynamic_cast<TH2F*> (fOutput->FindObject("fHistTrackTrigCorrDen"));
861 fHistTrigVtxCorrNum = dynamic_cast<TH2F*> (fOutput->FindObject("fHistTrigVtxCorrNum"));
862 fHistTrigVtxCorrDen = dynamic_cast<TH2F*> (fOutput->FindObject("fHistTrigVtxCorrDen"));
864 fHistTrigCorrDen = dynamic_cast<TH2F*> (fOutput->FindObject("fHistTrigCorrDen"));
867 fHistMCEtavsZTriggMCvtxEvts = dynamic_cast<TH2F*> (fOutput->FindObject("fHistMCEtavsZTriggMCvtxEvts"));
868 fHistMCEtavsZTriggESDvtxEvts = dynamic_cast<TH2F*> (fOutput->FindObject("fHistMCEtavsZTriggESDvtxEvts"));
869 fHistMCEtavsZ = dynamic_cast<TH2F*> (fOutput->FindObject("fHistMCEtavsZ"));
872 fHistTRradius = dynamic_cast<TH1F*> (fOutput->FindObject("fHistTRradius"));
874 fHistContributorsvsDeVtx = dynamic_cast<TH2F*> (fOutput->FindObject("fHistContributorsvsDeVtx"));
876 fHistoDetectableNotr = dynamic_cast<TH3F*> (fOutput->FindObject("fHistoDetectableNotr"));
877 fHistoDetectabletr = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectabletr"));
879 fHistoNonStoppingTracks = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoNonStoppingTracks"));
881 fHistoDetectedLay1 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectedLay1"));
882 fHistoDetectedLay2 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectedLay2"));
884 fHistoPt = dynamic_cast<TH1F*> (fOutput->FindObject("fHistoPt"));
886 fHistoDetectableTRm1 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTRm1"));
887 fHistoDetectableTR0 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR0"));
888 fHistoDetectableTR1 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR1"));
889 fHistoDetectableTR2 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR2"));
890 fHistoDetectableTR3 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR3"));
891 fHistoDetectableTR4 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR4"));
892 fHistoDetectableTR5 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR5"));
893 fHistoDetectableTR6 = dynamic_cast<TH2F*> (fOutput->FindObject("fHistoDetectableTR6"));
894 fHistoRTRm1 = dynamic_cast<TH1F*> (fOutput->FindObject("fHistoRTRm1"));