]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/EMCAL/AliEmcalClusterMaker.cxx
e0e5d43e369ef0daefa31295b58f7d2ffe449142
[u/mrichter/AliRoot.git] / PWG / EMCAL / AliEmcalClusterMaker.cxx
1 // $Id$
2 //
3 // Cluster maker task.
4 //
5 // Author: C.Loizides
6
7 #include <TChain.h>
8 #include <TClonesArray.h>
9 #include "AliAODCaloCluster.h"
10 #include "AliAODEvent.h"
11 #include "AliAnalysisManager.h"
12 #include "AliEMCALRecoUtils.h"
13 #include "AliESDCaloCluster.h"
14 #include "AliESDEvent.h"
15 #include "AliEmcalClusterMaker.h"
16
17 ClassImp(AliEmcalClusterMaker)
18
19 //________________________________________________________________________
20 AliEmcalClusterMaker::AliEmcalClusterMaker() : 
21   AliAnalysisTaskEmcal("AliEmcalClusterMaker", kFALSE),
22   fOutCaloName(),
23   fRecoUtils(0),
24   fEsdMode(kTRUE),
25   fOutClusters(0),
26   fEnergyDistBefore(0),
27   fEtaPhiDistBefore(0),
28   fEnergyTimeHistBefore(0),
29   fEnergyDistAfter(0),
30   fEtaPhiDistAfter(0),
31   fEnergyTimeHistAfter(0)
32 {
33   // Default constructor.
34 }
35
36 //________________________________________________________________________
37 AliEmcalClusterMaker::AliEmcalClusterMaker(const char *name, Bool_t histo) : 
38   AliAnalysisTaskEmcal(name, histo),
39   fOutCaloName("EmcClusters"),
40   fRecoUtils(0),
41   fEsdMode(kTRUE),
42   fOutClusters(0),
43   fEnergyDistBefore(0),
44   fEtaPhiDistBefore(0),
45   fEnergyTimeHistBefore(0),
46   fEnergyDistAfter(0),
47   fEtaPhiDistAfter(0),
48   fEnergyTimeHistAfter(0)
49 {
50   // Standard constructor.
51   
52   SetMakeGeneralHistograms(histo);
53
54   fBranchNames="ESD:AliESDRun.,AliESDHeader.,PrimaryVertex.";
55 }
56
57 //________________________________________________________________________
58 AliEmcalClusterMaker::~AliEmcalClusterMaker()
59 {
60   // Destructor
61 }
62
63 //________________________________________________________________________
64 void AliEmcalClusterMaker::UserCreateOutputObjects()
65 {
66   // Create my user objects.
67
68   AliAnalysisTaskEmcal::UserCreateOutputObjects();
69
70   if (fRecoUtils) {
71     fRecoUtils->InitNonLinearityParam();
72     fRecoUtils->Print("");
73   }
74
75   if (!fCreateHisto) return;
76
77   fEnergyDistBefore = new TH1F("hEnergyDistBefore","hEnergyDistBefore;E_{clus} (GeV)",1500,0,150);
78   fOutput->Add(fEnergyDistBefore);
79   fEtaPhiDistBefore = new TH2F("hEtaPhiDistBefore","hEtaPhiDistBefore;#eta;#phi",280,-0.7,0.7,800,1.3,3.3);
80   fOutput->Add(fEtaPhiDistBefore);
81   fEnergyTimeHistBefore = new TH2F("hEnergyTimeDistBefore","hEnergyTimeDistBefore;E_{clus} (GeV);time",60,0,30,500,0,1e-6);
82   fOutput->Add(fEnergyTimeHistBefore);
83   fEnergyDistAfter = new TH1F("hEnergyDistAfter","hEnergyDistAfter;E_{clus} (GeV)",1500,0,150);
84   fOutput->Add(fEnergyDistAfter);
85   fEtaPhiDistAfter = new TH2F("hEtaPhiDistAfter","hEtaPhiDistAfter;#eta;#phi",280,-0.7,0.7,800,1.3,3.3);
86   fOutput->Add(fEtaPhiDistAfter);
87   fEnergyTimeHistAfter = new TH2F("hEnergyTimeDistAfter","hEnergyTimeDistAfter;E_{clus} (GeV);time",60,0,30,500,0,1e-6);
88   fOutput->Add(fEnergyTimeHistAfter);
89   PostData(1, fOutput);
90 }
91
92 //________________________________________________________________________
93 void AliEmcalClusterMaker::ExecOnce() 
94 {
95   // Initialize the analysis.
96
97   // Do base class initializations and if it fails -> bail out
98   AliAnalysisTaskEmcal::ExecOnce();
99   if (!fInitialized) 
100     return;
101
102   if (dynamic_cast<AliAODEvent*>(InputEvent())) 
103     fEsdMode = kFALSE;
104
105   if (fEsdMode) 
106     fOutClusters = new TClonesArray("AliESDCaloCluster");
107   else 
108     fOutClusters = new TClonesArray("AliAODCaloCluster");
109
110   fOutClusters->SetName(fOutCaloName);
111
112   // post output in event if not yet present
113   if (!(InputEvent()->FindListObject(fOutCaloName))) {
114     InputEvent()->AddObject(fOutClusters);
115   } else {
116     fInitialized = kFALSE;
117     AliFatal(Form("%s: Container with same name %s already present. Aborting", GetName(), fOutCaloName.Data()));
118     return;
119   }
120 }
121
122 //________________________________________________________________________
123 Bool_t AliEmcalClusterMaker::Run() 
124 {
125   // Run the hadronic correction
126
127   // delete output
128   fOutClusters->Delete();
129
130   // loop over clusters
131   Int_t clusCount = 0;
132   Int_t entries   = fCaloClusters->GetEntries();
133   for (Int_t i=0; i<entries; ++i) {
134     AliVCluster *clus = static_cast<AliVCluster*>(fCaloClusters->At(i));
135     if (!clus || !clus->IsEMCAL())
136       continue;
137
138     if (fCreateHisto) {
139       fEnergyDistBefore->Fill(clus->E());
140       Float_t pos[3] ={0,0,0};
141       clus->GetPosition(pos);
142       TVector3 vec(pos);
143       fEtaPhiDistBefore->Fill(vec.Eta(),vec.Phi());
144       fEnergyTimeHistBefore->Fill(clus->E(),clus->GetTOF());
145     }
146
147     AliVCluster *oc = 0;
148     if (fEsdMode) {
149       AliESDCaloCluster *ec = dynamic_cast<AliESDCaloCluster*>(clus);
150       if (!ec) continue;
151       oc = new ((*fOutClusters)[clusCount]) AliESDCaloCluster(*ec);
152     } else { 
153       AliAODCaloCluster *ac = dynamic_cast<AliAODCaloCluster*>(clus);
154       if (!ac) continue;
155       oc = new ((*fOutClusters)[clusCount]) AliAODCaloCluster(*ac);
156     }
157     if (fRecoUtils) {
158       if (fRecoUtils->IsRejectExoticCluster()) {
159         if (fRecoUtils->IsExoticCluster(oc,fCaloCells))
160         continue;
161       }
162       if (fRecoUtils->GetNonLinearityFunction()!=AliEMCALRecoUtils::kNoCorrection) {
163         Double_t energy = fRecoUtils->CorrectClusterEnergyLinearity(oc);
164         oc->SetE(energy);
165       }
166     }
167     if (!AcceptCluster(oc))
168       continue;
169     clusCount++;
170
171     if (fCreateHisto) {
172       fEnergyDistAfter->Fill(oc->E());
173       Float_t pos[3] ={0,0,0};
174       oc->GetPosition(pos);
175       TVector3 vec(pos);
176       fEtaPhiDistAfter->Fill(vec.Eta(),vec.Phi());
177       fEnergyTimeHistAfter->Fill(oc->E(),oc->GetTOF());
178     }
179
180   }
181   if ((clusCount>0) && (clusCount==fOutClusters->GetEntries()))
182     fOutClusters->RemoveAt(clusCount);
183   return kTRUE;
184 }