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