]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetHistos.cxx
cancelling accidental remove of 'delete anaUtil'
[u/mrichter/AliRoot.git] / JETAN / AliJetHistos.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
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 /* $Id$ */
17
18 //---------------------------------------------------------------------
19 // Jet Histos class
20 // Creates and fills a few cummon histograms for jet analysis
21 //
22 //---------------------------------------------------------------------
23
24 #include <TList.h>
25 #include <TClonesArray.h>
26 #include <TH1I.h>
27 #include <TH1F.h>
28 #include <TMath.h>
29
30 #include "AliAODJet.h"
31 #include "AliJetHistos.h"
32
33 ClassImp(AliJetHistos)
34
35 ///////////////////////////////////////////////////////////////////////
36
37 AliJetHistos::AliJetHistos():
38   fNJetsH(0x0),
39   fPtH(0x0),
40   fEtaH(0x0),
41   fEneH(0x0),
42   fPhiH(0x0)
43 {
44   // Default constructor
45 }
46
47 //-----------------------------------------------------------------------
48 void AliJetHistos::CreateHistos()
49 {
50   // create histos
51
52   fNJetsH = new TH1I("NJetsH","Number of Jets",12,0,11);
53   SetProperties(fNJetsH,"Number of jets","Entries");
54
55   fPtH = new TH1F("PtH","Pt of Jets",50,0.,200.);
56   SetProperties(fPtH,"P_{#perp} [GeV]","Entries");
57
58   fEtaH = new TH1F("EtaH","Pseudorapidity of Jets",30,-1.5,1.5);
59   SetProperties(fEtaH,"#eta","Entries");
60
61   fEneH = new TH1F("EneH","Energy of Jets",50,0.,200.);
62   SetProperties(fEneH,"Energy [GeV]","Entries");
63
64   fPhiH = new TH1F("PhiH","Azimuthal angle of Jets",
65                    60,0.,2.0*TMath::Pi());
66   SetProperties(fPhiH,"#phi","Entries");
67 }
68
69 //-----------------------------------------------------------------------
70 AliJetHistos::~AliJetHistos()
71 {
72   // Destructor
73   delete fNJetsH;
74   delete fPtH;
75   delete fEtaH;
76   delete fEneH;
77   delete fPhiH;
78 }
79
80 //-----------------------------------------------------------------------
81 void AliJetHistos::SetProperties(TH1* h,const char* x, const char* y) const
82 {
83   // Sets the histogram style properties
84   h->SetMarkerStyle(20);
85   h->SetMarkerSize(.5);
86   h->SetMarkerColor(2);
87   h->SetXTitle(x);
88   h->SetYTitle(y);
89   h->Sumw2();
90 }
91
92 //-----------------------------------------------------------------------
93 void AliJetHistos::AddHistosToList(TList *list) const
94 {
95   // Add histos to the list
96   list->Add(fNJetsH);
97   list->Add(fPtH);
98   list->Add(fEtaH);
99   list->Add(fEneH);
100   list->Add(fPhiH);
101 }
102
103 //-----------------------------------------------------------------------
104 void AliJetHistos::FillHistos(TClonesArray *jets)
105 {
106   // Fill histograms
107   if(!jets)return;
108   Int_t nj = jets->GetEntries();
109   fNJetsH->Fill(nj,1);
110
111   if (nj == 0 ) return;
112
113   AliAODJet *j;
114   for (Int_t i=0;i<nj;i++) {
115     j = (AliAODJet *) jets->At(i);
116     fPtH->Fill(j->Pt(),1);
117     fEtaH->Fill(j->Eta(),1);
118     fEneH->Fill(j->E(),1);
119     fPhiH->Fill(j->Phi(),1);
120   }
121   
122 }