]> git.uio.no Git - u/mrichter/AliRoot.git/blame - JETAN/AliJetControlPlots.cxx
Fast EMCAL simulation option added.
[u/mrichter/AliRoot.git] / JETAN / AliJetControlPlots.cxx
CommitLineData
99e5fe42 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//---------------------------------------------------------------------
17// Jet Control Plots class
18// manages histograms with control plots of jet searching
19// Stores the output of a jet algorithm
20// Author: jgcn@mda.cinvestav.mx
21//---------------------------------------------------------------------
22
23#include <TStyle.h>
24#include <TCanvas.h>
25#include <TLorentzVector.h>
26#include <TFile.h>
27#include <TClonesArray.h>
28#include <TH1I.h>
29#include <TH1D.h>
30#include "AliJetReader.h"
31#include "AliJet.h"
32
33#include "AliJetControlPlots.h"
34ClassImp(AliJetControlPlots)
35
36////////////////////////////////////////////////////////////////////////
37
38AliJetControlPlots::AliJetControlPlots()
39{
40 //
41 // Constructor
42 //
43 fNJetT=0;
44
45 fNJetsH = new TH1I("fNJetsH","Number of Jets",12,0,11);
46 SetProperties(fNJetsH,"Number of jets","Entries");
47
48 fMultH = new TH1I("fMultH","Multiplicity of Jets",22,0,21);
49 SetProperties(fMultH,"Multiplicity of jets","Entries");
50
51 fPtH = new TH1D("fPtH","Pt of Jets",50,0.,200.);
52 SetProperties(fPtH,"P_{t} [GeV]","Entries");
53
54 fEtaH = new TH1D("fEtaH","Pseudorapidity of Jets",30,-1.5,1.5);
55 SetProperties(fEtaH,"#eta","Entries");
56
57 fEneH = new TH1D("fEneH","Energy of Jets",50,0.,200.);
58 SetProperties(fEneH,"Energy [GeV]","Entries");
59
60 fFragH = new TH1D("fFragH","Jet Fragmentation",20,0.,1.);
61 SetProperties(fFragH,"x=E_{i}/E_{JET}","1/N_{JET}dN_{ch}/dx");
62
63 fFragLnH = new TH1D("fFragLnH","Jet Fragmentation",20,0.,10);
64 SetProperties(fFragLnH,"#xi=ln(E_{JET}/E_{i}","1/N_{JET}dN_{ch}/d#xi");
65
66 fPhiH = new TH1D("fPhiH","Azimuthal angle of Jets",
67 60,-TMath::Pi(),TMath::Pi());
68 SetProperties(fPhiH,"#phi","Entries");
69
70 fInJetH = new TH1D("fInJetH","Percentage of particles in jets",
71 50,0,1);
72 SetProperties(fInJetH,"percentage of particles in jets","Entries");
73}
74
75////////////////////////////////////////////////////////////////////////
76
77AliJetControlPlots::~AliJetControlPlots()
78{
79 //
80 // Destructor
81 //
82 delete fNJetsH;
83 delete fMultH;
84 delete fPtH;
85 delete fEtaH;
86 delete fEneH;
87 delete fFragH;
88 delete fFragLnH;
89 delete fPhiH;
90 delete fInJetH;
91}
92
93////////////////////////////////////////////////////////////////////////
94
95void AliJetControlPlots::FillHistos(AliJet *j, AliJetReader *r)
96{
97// Fills the histograms
98
99 Int_t nj = j->GetNJets();
100 fNJetsH->Fill(nj);
101 fNJetT+=nj;
102
103 // kinematics, occupancy and multiplicities
104 TArrayI mj = j->GetMultiplicities();
105 Int_t mjTot=0;
106 for (Int_t i=0;i<nj;i++) {
107 mjTot+=mj[i];
108 fMultH->Fill(mj[i]);
109 fPtH->Fill(j->GetPt(i));
110 fEtaH->Fill(j->GetEta(i));
111 fEneH->Fill(j->GetE(i));
112 fPhiH->Fill(j->GetPhi(i));
113 }
114 if (nj>0)
115 fInJetH->Fill(((Double_t) mjTot)/((Double_t) j->GetNinput()));
116
117 // fragmentation
118 TClonesArray *lvArray = r->GetMomentumArray();
119 TArrayI inJet = j->GetInJet();
120 Int_t nIn = j->GetNinput();
121 for(Int_t i=0;i<nIn;i++) {
122 if (inJet[i] == -1) continue;
123 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
124 Double_t xe = (lv->E())/(j->GetE(inJet[i]));
125 fFragH->Fill(xe);
126 fFragLnH->Fill(TMath::Log(1.0/xe));
127 }
128}
129
130////////////////////////////////////////////////////////////////////////
131
132void AliJetControlPlots::Normalize()
133{
134 if (fNJetT == 0) return;
135 fFragH->Sumw2();
136 fFragH->Scale(20.0/((Double_t) fNJetT));
137 fFragLnH->Sumw2();
138 fFragLnH->Scale(2.0/((Double_t) fNJetT));
139}
140
141////////////////////////////////////////////////////////////////////////
142
143void AliJetControlPlots::PlotHistos()
144{
145 // style
146 gStyle->SetOptStat(kFALSE);
147 gStyle->SetOptTitle(kFALSE);
148
149 TCanvas *c = new TCanvas("c","AliJetControlPlots",50,200,900,700);
150 c->Divide(3,3);
151 c->cd(1); fNJetsH->Draw("e1");
152 c->cd(2); fMultH->Draw("e1");
153 c->cd(3); fInJetH->Draw("e1");
154 c->cd(4); fPtH->Draw("e1");
155 c->cd(5); fEtaH->Draw("e1");
156 c->cd(6); fPhiH->Draw("e1");
157 c->cd(7); fEneH->Draw("e1");
158 c->cd(8); fFragH->Draw("e1");
159 c->cd(9); fFragLnH->Draw("e1");
160}
161
162////////////////////////////////////////////////////////////////////////
163
164void AliJetControlPlots::SetProperties(TH1* h,const char* x, const char* y) const
165{
166//
167// Sets the histogram style properties
168 h->SetMarkerStyle(20);
169 h->SetMarkerSize(.5);
170 h->SetMarkerColor(2);
171 h->SetXTitle(x);
172 h->SetYTitle(y);
173}
174
175