1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 //---------------------------------------------------------------------
27 #include "AliJetReader.h"
30 #include "AliJetControlPlots.h"
31 ClassImp(AliJetControlPlots)
33 ////////////////////////////////////////////////////////////////////////
35 AliJetControlPlots::AliJetControlPlots()
41 fNJetsH = new TH1I("fNJetsH","Number of Jets",12,0,11);
42 SetProperties(fNJetsH,"Number of jets","Entries");
44 fMultH = new TH1I("fMultH","Multiplicity of Jets",22,0,21);
45 SetProperties(fMultH,"Multiplicity of jets","Entries");
47 fInJetH = new TH1D("fInJetH","Percentage of particles in jets",50,0,1);
48 SetProperties(fInJetH,"percentage of particles in jets","Entries");
51 fPtH = new TH1D("fPtH","Pt of Jets",50,0.,200.);
52 SetProperties(fPtH,"P_{t} [GeV]","Entries");
54 fEtaH = new TH1D("fEtaH","Pseudorapidity of Jets",30,-1.5,1.5);
55 SetProperties(fEtaH,"#eta","Entries");
57 fEneH = new TH1D("fEneH","Energy of Jets",50,0.,200.);
58 SetProperties(fEneH,"Energy [GeV]","Entries");
60 fPhiH = new TH1D("fPhiH","Azimuthal angle of Jets",
61 60,0.,2.0*TMath::Pi());
62 SetProperties(fPhiH,"#phi","Entries");
65 fFragH = new TH1D("fFragH","Leading Jet Fragmentation (selected tracks)",
67 SetProperties(fFragH,"x=E_{i}/E_{JET}","1/N_{JET}dN_{ch}/dx");
69 fFragrH = new TH1D("fFragrH","Leading Jet Fragmentation (rejected tracks)",
71 SetProperties(fFragrH,"x=E_{i}/E_{JET}","1/N_{JET}dN_{ch}/dx");
74 fFragLnH = new TH1D("fFragLnH","Leading Jet Fragmentation (selected tracks)",
76 SetProperties(fFragLnH,"#xi=ln(E_{JET}/E_{i})","1/N_{JET}dN_{ch}/d#xi");
78 fFragLnrH = new TH1D("fFragLnrH",
79 "Leading Jet Fragmentation (rejected tracks)",
81 SetProperties(fFragLnrH,"#xi=ln(E_{JET}/E_{i})","1/N_{JET}dN_{ch}/d#xi");
84 fShapeH = new TH1D("fShapeH","Leading Jet Shape (selected tracks)",
86 SetProperties(fShapeH,"r","1/N_{JET}#Sigma P_{T}(0,r)/P_{T}(0,R)");
88 fShaperH = new TH1D("fShaperH","Leading Jet Shape (rejected tracks)",
90 SetProperties(fShaperH,"r","1/N_{JET}#Sigma P_{T}(0,r)/P_{T}(0,R)");
93 ////////////////////////////////////////////////////////////////////////
95 AliJetControlPlots::~AliJetControlPlots()
113 ////////////////////////////////////////////////////////////////////////
115 void AliJetControlPlots::FillHistos(AliJet *j)
117 // Fills the histograms
118 Int_t nj = j->GetNJets();
122 // kinematics, occupancy and multiplicities
123 TArrayI mj = j->GetMultiplicities();
125 for (Int_t i=0;i<nj;i++) {
127 fMultH->Fill(mj[i],1);
128 fPtH->Fill(j->GetPt(i),1);
129 fEtaH->Fill(j->GetEta(i),1);
130 fEneH->Fill(j->GetE(i),1);
131 fPhiH->Fill(j->GetPhi(i),1);
133 fInJetH->Fill(((Double_t) mjTot)/((Double_t) j->GetNinput()),1);
135 // fragmentation of leading jet
136 TArrayI inJet = j->GetInJet();
137 TArrayF etain = j->GetEtaIn();
138 TArrayF ptin = j->GetPtIn();
139 for(Int_t i=0;i<(inJet.GetSize());i++) {
140 Float_t t1 = TMath::Tan(2.0*TMath::ATan(TMath::Exp(-etain[i])));
141 Float_t ene = ptin[i]*TMath::Sqrt(1.+1./(t1*t1));
143 fFragH->Fill((Float_t) ene/(j->GetE(0)),1);
144 fFragLnH->Fill((Float_t) TMath::Log((j->GetE(0))/ene),1);
146 if (inJet[i] == -1) {
147 fFragrH->Fill((Float_t) ene/(j->GetE(0)),1);
148 fFragLnrH->Fill((Float_t) TMath::Log((j->GetE(0))/ene),1);
152 // shape of leading jet
153 // * CAREFUL: depends on number of bins and bin sizes
154 // of shape histo. HARDWIRED at the moment *
155 // Plot against dr NOT dr/R_jet!!!
156 TArrayF phiin = j->GetPhiIn();
157 Float_t rptS[20], rptR[20];
158 for(Int_t i=0;i<20;i++) rptS[i]=rptR[i]=0.0;
160 for(Int_t i=0;i<inJet.GetSize();i++) {
161 if (inJet[i] == 1 || inJet[i] == -1) {
162 Float_t deta = etain[i] - j->GetEta(0);
163 Float_t dphi = phiin[i] - j->GetPhi(0);
164 if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
165 if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
166 Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
168 Int_t bin = (Int_t) TMath::Floor(dr/0.05);
169 if (inJet[i] == 1) rptS[bin]+=ptin[i]/(j->GetPt(0));
170 if (inJet[i] == -1) rptR[bin]+=ptin[i]/(j->GetPt(0));
176 for (Int_t i=0;i<20;i++) {
180 fShapeH->Fill(r,ptS);
181 fShaperH->Fill(r,ptR);
186 ////////////////////////////////////////////////////////////////////////
188 void AliJetControlPlots::Normalize()
190 // *CAREFUL: depends on histogram number of bins
191 if (fNJetT == 0) return;
192 fFragH->Scale(20.0/((Double_t) fNJetT));
193 fFragLnH->Scale(2.0/((Double_t) fNJetT));
194 fFragrH->Scale(20.0/((Double_t) fNJetT));
195 fFragLnrH->Scale(2.0/((Double_t) fNJetT));
196 fShapeH->Scale(1.0/((Double_t) fNJetT));
197 fShaperH->Scale(1.0/((Double_t) fNJetT));
200 ////////////////////////////////////////////////////////////////////////
202 void AliJetControlPlots::PlotHistos()
205 gStyle->SetOptStat(kFALSE);
206 gStyle->SetOptTitle(kFALSE);
208 TCanvas *c = new TCanvas("c","AliJetControlPlots",50,200,1200,700);
210 c->cd(1); fNJetsH->Draw("e1");
211 c->cd(2); fMultH->Draw("e1");
212 c->cd(3); fInJetH->Draw("e1");
213 c->cd(4); fPtH->Draw("e1");
214 c->cd(5); fEtaH->Draw("e1");
215 c->cd(6); fPhiH->Draw("e1");
216 c->cd(7); fEneH->Draw("e1");
217 c->cd(8); fFragLnH->Draw("e1");
218 c->cd(9); fFragLnrH->Draw("e1");
219 c->cd(10); fShapeH->Draw("e1");
220 c->cd(11); fShaperH->Draw("e1");
223 ////////////////////////////////////////////////////////////////////////
225 void AliJetControlPlots::SetProperties(TH1* h,const char* x, const char* y) const
227 // Sets the histogram style properties
228 h->SetMarkerStyle(20);
229 h->SetMarkerSize(.5);
230 h->SetMarkerColor(2);