]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/physics/AliHLTCaloHistoClusterEnergy.cxx
Add output of VZERO reco directly to multCorr
[u/mrichter/AliRoot.git] / HLT / global / physics / AliHLTCaloHistoClusterEnergy.cxx
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        * 
3  * All rights reserved.                                                   *
4  *                                                                        *
5  * Primary Authors: Svein Lindal <slindal@fys.uio.no>                     *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          * 
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /** 
17  * @file   AliHLTCaloHistoClusterEnergy
18  * @author Svein Lindal
19  * @date 
20  * @brief  Produces histograms of cluster energy distributions 
21  */
22
23 // see header file for class documentation
24 // or
25 // refer to README to build package
26 // or
27 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
28
29 #include "AliHLTCaloHistoClusterEnergy.h"
30 #include "AliESDCaloCluster.h"
31 #include "AliHLTCaloClusterDataStruct.h"
32 #include "TH1F.h"
33 #include "TH2F.h"
34 #include "TRefArray.h"
35 #include "TVector3.h"
36 #include "TMath.h"
37
38 ClassImp(AliHLTCaloHistoClusterEnergy);
39
40 AliHLTCaloHistoClusterEnergy::AliHLTCaloHistoClusterEnergy(TString det) :
41   AliHLTCaloHistoProducer(),
42   fHistClusterEnergy(NULL),
43   fHistClusterEnergyVsNCells(NULL),
44   fHistClusterEnergyDepositEtaPhi(NULL)
45 {
46   // See header file for documentation
47   fHistClusterEnergy = new TH1F(Form("%s fHistClusterEnergy", det.Data()), Form("%s Distribution of total energy in clusters", det.Data()), 5000, 0, 100);
48   fHistClusterEnergy->GetXaxis()->SetTitle("E GeV");
49   fHistClusterEnergy->GetYaxis()->SetTitle("Number of counts");
50   fHistArray->AddLast(fHistClusterEnergy);
51
52   fHistClusterEnergyVsNCells = new TH2F(Form("%s fHistClusterEnergyVsNCells", det.Data()), Form("%s Distribution of Energy vs Number of Cells in cluster", det.Data()), 1000, 0, 100, 50, 0 , 50);
53   fHistClusterEnergyVsNCells->GetXaxis()->SetTitle("Energy in cluster (GeV)");
54   fHistClusterEnergyVsNCells->GetYaxis()->SetTitle("Number of Cells in cluster");
55   fHistArray->AddLast(fHistClusterEnergyVsNCells);
56   
57   Float_t phiMin = 0.;
58   Float_t phiMax = 2*TMath::Pi();
59   Float_t etaMin = -1.;
60   Float_t etaMax = 1.;
61   
62   if(det == "PHOS")
63   {
64      phiMin = 255.0/180.*TMath::Pi();
65      phiMax = 325.0/180.*TMath::Pi();
66      etaMin = -0.13;
67      etaMax = 0.13;
68   }
69   
70   fHistClusterEnergyDepositEtaPhi = new TH2F(Form("%s fHistClusterEnergyDepositedEtaPhi", det.Data()), Form("%s Amount of energy deposited in Phi vs Eta", det.Data()), 200, phiMin, phiMax, 50, etaMin , etaMax);
71   fHistClusterEnergyDepositEtaPhi->GetXaxis()->SetTitle("#phi");
72   fHistClusterEnergyDepositEtaPhi->GetYaxis()->SetTitle("#eta");
73   fHistArray->AddLast(fHistClusterEnergyDepositEtaPhi);
74
75 }
76
77 AliHLTCaloHistoClusterEnergy::~AliHLTCaloHistoClusterEnergy()
78 {
79   //destructor
80   if(fHistClusterEnergy)
81     delete fHistClusterEnergy;
82   fHistClusterEnergy = NULL;
83
84   if(fHistClusterEnergyVsNCells)
85     delete fHistClusterEnergyVsNCells;
86   fHistClusterEnergyVsNCells = NULL;
87 }
88
89 Int_t AliHLTCaloHistoClusterEnergy::FillHistograms(Int_t nc, TRefArray * clusterArray) {
90   //See header file for documentation
91
92   for(int ic = 0; ic < nc; ic++) {
93     AliESDCaloCluster * cluster = static_cast<AliESDCaloCluster*>(clusterArray->At(ic));
94     FillClusterEnergyHistos(cluster);
95   }
96   return 0;
97 }
98
99 Int_t AliHLTCaloHistoClusterEnergy::FillHistograms(Int_t nc, vector<AliHLTCaloClusterDataStruct*> &cVec) {
100   //See header file for documentation
101   // HLTInfo("histo");
102   for(int ic = 0; ic < nc; ic++) {
103     AliHLTCaloClusterDataStruct * cluster = cVec.at(ic);
104     FillClusterEnergyHistos(cluster);
105   }
106   return 0;
107 }
108
109 template <class T>
110 Int_t AliHLTCaloHistoClusterEnergy::FillClusterEnergyHistos(T* cluster) {
111   fHistClusterEnergy->Fill(cluster->E());
112   fHistClusterEnergyVsNCells->Fill(cluster->E(), cluster->GetNCells());
113   
114   Float_t pos[3];
115   cluster->GetPosition(pos);
116   TVector3 vec(pos);
117   
118   // Stupid hack, too tired to fix
119   if(vec.Phi() < 0)  fHistClusterEnergyDepositEtaPhi->Fill(2*TMath::Pi() + vec.Phi(), vec.Eta(), cluster->E());
120   else fHistClusterEnergyDepositEtaPhi->Fill(vec.Phi(), vec.Eta(), cluster->E());
121   
122   return 0;
123 }