]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/totEt/AliAnalysisTaskTotEt.cxx
- adding cut on probability of PID for histograms concerning identified particles
[u/mrichter/AliRoot.git] / PWG4 / totEt / AliAnalysisTaskTotEt.cxx
CommitLineData
cf6522d1 1//_________________________________________________________________________
2// Utility Class for transverse energy studies
3// Task for analysis
4// - reconstruction and MC output
5// implementation file
6//
7//*-- Authors: Oystein Djuvsland (Bergen), David Silvermyr (ORNL)
8//_________________________________________________________________________
9
2fbf38ac 10#include "TChain.h"
cf6522d1 11#include "TList.h"
2fbf38ac 12#include "TH2F.h"
2fbf38ac 13
14#include "AliESDEvent.h"
2fbf38ac 15#include "AliMCEvent.h"
cf6522d1 16
2fbf38ac 17#include "AliAnalysisTaskTotEt.h"
18#include "AliAnalysisEtReconstructedPhos.h"
19#include "AliAnalysisEtReconstructedEmcal.h"
20#include "AliAnalysisEtMonteCarloPhos.h"
21#include "AliAnalysisEtMonteCarloEmcal.h"
22
23#include <iostream>
2fbf38ac 24
25using namespace std;
26
27ClassImp(AliAnalysisTaskTotEt)
28
29//________________________________________________________________________
30AliAnalysisTaskTotEt::AliAnalysisTaskTotEt(const char *name) :
31 AliAnalysisTaskSE(name)
2fbf38ac 32 ,fOutputList(0)
33 ,fRecAnalysis(0)
34 ,fMCAnalysis(0)
35 ,fHistEtRecvsEtMC(0)
2fbf38ac 36{
37 // Constructor
38
39 // select if we should use EMCal or PHOS class
40 // PHOS by default, EMCal if name string contains EMC
41 TString t(name);
42 t.ToUpper();
43 if (t.Contains("EMC")) {
44 fRecAnalysis = new AliAnalysisEtReconstructedEmcal();
45 fMCAnalysis = new AliAnalysisEtMonteCarloEmcal();
46 }
47 else {
48 fRecAnalysis = new AliAnalysisEtReconstructedPhos();
49 fMCAnalysis = new AliAnalysisEtMonteCarloPhos();
50 }
51
52 fRecAnalysis->Init();
53 fMCAnalysis->Init();
54
2fbf38ac 55 // Define input and output slots here
56 // Input slot #0 works with a TChain
57 DefineInput(0, TChain::Class());
58 // Output slot #1 writes into a TH1 container
59
60 DefineOutput(1, TList::Class());
61
62}
63
64
65//________________________________________________________________________
66void AliAnalysisTaskTotEt::UserCreateOutputObjects()
67{
68 // Create histograms
69 // Called once
70 fMCAnalysis->CreateHistograms();
71 fRecAnalysis->CreateHistograms();
72 fOutputList = new TList;
73 fRecAnalysis->FillOutputList(fOutputList);
74 fMCAnalysis->FillOutputList(fOutputList);
75 fHistEtRecvsEtMC = new TH2F("fHistEtRecvsEtMC", "Reconstructed E_{t} vs MC E_{t}", 1000, 0.000, 100, 1000, 0.0001, 100);
76 fOutputList->Add(fHistEtRecvsEtMC);
77}
78
79//________________________________________________________________________
80void AliAnalysisTaskTotEt::UserExec(Option_t *)
cf6522d1 81{ // execute method
82 AliESDEvent *event = dynamic_cast<AliESDEvent*>(InputEvent());
2fbf38ac 83 if (!event) {
84 Printf("ERROR: Could not retrieve event");
85 return;
86 }
87
88 fRecAnalysis->AnalyseEvent(event);
89
90 AliMCEvent* mcEvent = MCEvent();
91 if (mcEvent)
92 {
93 fMCAnalysis->AnalyseEvent(mcEvent);
94 }
95
96 fHistEtRecvsEtMC->Fill(fRecAnalysis->GetTotEtAcc(), fMCAnalysis->GetTotEt());
97
98// Post output data.
99 PostData(1, fOutputList);
100
101}
102
103//________________________________________________________________________
104void AliAnalysisTaskTotEt::Terminate(Option_t *)
105{
106 // Draw result to the screen
107 // Called once at the end of the query
108
109 fOutputList = dynamic_cast<TList*> (GetOutputData(1));
110 if (!fOutputList) {
111 printf("ERROR: Output list not available\n");
112 return;
113 }
114}
115
116
117