]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/totEt/AliAnalysisTaskTotEt.cxx
new classes for tracking efficiency
[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"
b5821c13 16#include "AliESDtrackCuts.h"
cf6522d1 17
2fbf38ac 18#include "AliAnalysisTaskTotEt.h"
19#include "AliAnalysisEtReconstructedPhos.h"
20#include "AliAnalysisEtReconstructedEmcal.h"
21#include "AliAnalysisEtMonteCarloPhos.h"
22#include "AliAnalysisEtMonteCarloEmcal.h"
23
24#include <iostream>
2fbf38ac 25
26using namespace std;
27
28ClassImp(AliAnalysisTaskTotEt)
29
30//________________________________________________________________________
31AliAnalysisTaskTotEt::AliAnalysisTaskTotEt(const char *name) :
32 AliAnalysisTaskSE(name)
2fbf38ac 33 ,fOutputList(0)
34 ,fRecAnalysis(0)
35 ,fMCAnalysis(0)
36 ,fHistEtRecvsEtMC(0)
b5821c13 37 ,fEsdtrackCutsTPC(0)
2fbf38ac 38{
39 // Constructor
40
41 // select if we should use EMCal or PHOS class
42 // PHOS by default, EMCal if name string contains EMC
43 TString t(name);
44 t.ToUpper();
45 if (t.Contains("EMC")) {
46 fRecAnalysis = new AliAnalysisEtReconstructedEmcal();
47 fMCAnalysis = new AliAnalysisEtMonteCarloEmcal();
48 }
49 else {
50 fRecAnalysis = new AliAnalysisEtReconstructedPhos();
51 fMCAnalysis = new AliAnalysisEtMonteCarloPhos();
52 }
53
54 fRecAnalysis->Init();
55 fMCAnalysis->Init();
56
2fbf38ac 57 // Define input and output slots here
58 // Input slot #0 works with a TChain
59 DefineInput(0, TChain::Class());
60 // Output slot #1 writes into a TH1 container
61
62 DefineOutput(1, TList::Class());
63
64}
464aa50c 65AliAnalysisTaskTotEt::~AliAnalysisTaskTotEt(){//Destructor
951efd81 66 fOutputList->Clear();
67 delete fOutputList;
68 delete fRecAnalysis;
69 delete fMCAnalysis;
70 delete fEsdtrackCutsTPC;
71}
2fbf38ac 72
73//________________________________________________________________________
74void AliAnalysisTaskTotEt::UserCreateOutputObjects()
75{
76 // Create histograms
77 // Called once
78 fMCAnalysis->CreateHistograms();
79 fRecAnalysis->CreateHistograms();
80 fOutputList = new TList;
951efd81 81 fOutputList->SetOwner();
2fbf38ac 82 fRecAnalysis->FillOutputList(fOutputList);
83 fMCAnalysis->FillOutputList(fOutputList);
84 fHistEtRecvsEtMC = new TH2F("fHistEtRecvsEtMC", "Reconstructed E_{t} vs MC E_{t}", 1000, 0.000, 100, 1000, 0.0001, 100);
85 fOutputList->Add(fHistEtRecvsEtMC);
b5821c13 86
87
88 fEsdtrackCutsTPC = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
89 fEsdtrackCutsTPC->SetName("fEsdTrackCutsTPCOnly");
90 fOutputList->Add(fEsdtrackCutsTPC);
91 if(fEsdtrackCutsTPC){
92 fRecAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
93 fMCAnalysis->SetTPCOnlyTrackCuts( GetTPCOnlyTrackCuts());
94 }
95 else{
96 Printf("Error: no track cuts!");
97 }
98
2fbf38ac 99}
100
101//________________________________________________________________________
102void AliAnalysisTaskTotEt::UserExec(Option_t *)
cf6522d1 103{ // execute method
104 AliESDEvent *event = dynamic_cast<AliESDEvent*>(InputEvent());
2fbf38ac 105 if (!event) {
106 Printf("ERROR: Could not retrieve event");
107 return;
108 }
109
110 fRecAnalysis->AnalyseEvent(event);
111
112 AliMCEvent* mcEvent = MCEvent();
113 if (mcEvent)
114 {
115 fMCAnalysis->AnalyseEvent(mcEvent);
116 }
117
118 fHistEtRecvsEtMC->Fill(fRecAnalysis->GetTotEtAcc(), fMCAnalysis->GetTotEt());
119
120// Post output data.
121 PostData(1, fOutputList);
122
123}
124
125//________________________________________________________________________
126void AliAnalysisTaskTotEt::Terminate(Option_t *)
127{
128 // Draw result to the screen
129 // Called once at the end of the query
130
131 fOutputList = dynamic_cast<TList*> (GetOutputData(1));
132 if (!fOutputList) {
133 printf("ERROR: Output list not available\n");
134 return;
135 }
136}
137
138
139