]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetControlPlots.cxx
Protection added.
[u/mrichter/AliRoot.git] / JETAN / AliJetControlPlots.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 Control Plots class 
20 // manages histograms with control plots of jet searching
21 // Stores the output of a jet algorithm
22 // Author: jgcn@mda.cinvestav.mx
23 //---------------------------------------------------------------------
24
25 #include <TCanvas.h>
26 #include <TH1D.h>
27 #include <TH1I.h>
28 #include <TMath.h>
29 #include <TStyle.h>
30
31 #include "AliJetReader.h"
32 #include "AliJet.h"
33
34 #include "AliJetControlPlots.h"
35 ClassImp(AliJetControlPlots)
36   
37 ////////////////////////////////////////////////////////////////////////
38
39 AliJetControlPlots::AliJetControlPlots():
40   fNJetsH(0),
41   fMultH(0),
42   fPtH(0),
43   fEtaH(0),
44   fEneH(0),
45   fFragH(0),
46   fFragLnH(0),
47   fFragrH(0),
48   fFragLnrH(0),
49   fShapeH(0),
50   fShaperH(0),
51   fPhiH(0),
52   fInJetH(0),
53   fNJetT(0)
54
55 {
56   // Default constructor
57
58   // general properties
59   fNJetsH = new TH1I("fNJetsH","Number of Jets",12,0,11);
60   SetProperties(fNJetsH,"Number of jets","Entries");
61
62   fMultH = new TH1I("fMultH","Multiplicity of Jets",22,0,21);
63   SetProperties(fMultH,"Multiplicity of jets","Entries");
64
65   fInJetH = new  TH1D("fInJetH","Percentage of particles in jets",50,0,1);
66   SetProperties(fInJetH,"percentage of particles in jets","Entries");
67
68   // kinematics
69   fPtH = new TH1D("fPtH","Pt of Jets",50,0.,200.);
70   SetProperties(fPtH,"P_{t} [GeV]","Entries");
71
72   fEtaH = new TH1D("fEtaH","Pseudorapidity of Jets",30,-1.5,1.5);
73   SetProperties(fEtaH,"#eta","Entries");
74
75   fEneH = new TH1D("fEneH","Energy of Jets",50,0.,200.);
76   SetProperties(fEneH,"Energy [GeV]","Entries");
77
78   fPhiH = new TH1D("fPhiH","Azimuthal angle of Jets",
79                    60,0.,2.0*TMath::Pi());
80   SetProperties(fPhiH,"#phi","Entries");
81
82   // fragmentation 
83   fFragH = new TH1D("fFragH","Leading Jet Fragmentation (selected tracks)",
84                     20,0.,1.);
85   SetProperties(fFragH,"x=E_{i}/E_{JET}","1/N_{JET}dN_{ch}/dx");
86
87   fFragrH = new TH1D("fFragrH","Leading Jet Fragmentation (rejected tracks)",
88                     20,0.,1.);
89   SetProperties(fFragrH,"x=E_{i}/E_{JET}","1/N_{JET}dN_{ch}/dx");
90
91   // fragmentation log
92   fFragLnH = new TH1D("fFragLnH","Leading Jet Fragmentation (selected tracks)",
93                       20,0.,10);
94   SetProperties(fFragLnH,"#xi=ln(E_{JET}/E_{i})","1/N_{JET}dN_{ch}/d#xi");
95
96   fFragLnrH = new TH1D("fFragLnrH",
97                        "Leading Jet Fragmentation (rejected tracks)",
98                        20,0.,10);
99   SetProperties(fFragLnrH,"#xi=ln(E_{JET}/E_{i})","1/N_{JET}dN_{ch}/d#xi");
100
101   // jet shape
102   fShapeH = new TH1D("fShapeH","Leading Jet Shape (selected tracks)",
103                      20,0.,1.);
104   SetProperties(fShapeH,"r","1/N_{JET}#Sigma P_{T}(0,r)/P_{T}(0,R)");
105       
106   fShaperH = new TH1D("fShaperH","Leading Jet Shape (rejected tracks)",
107                      20,0.,1.);
108   SetProperties(fShaperH,"r","1/N_{JET}#Sigma P_{T}(0,r)/P_{T}(0,R)");
109 }
110
111 ////////////////////////////////////////////////////////////////////////
112
113 AliJetControlPlots::~AliJetControlPlots()
114 {
115   // Destructor
116   delete fNJetsH;
117   delete fMultH;
118   delete fPtH;
119   delete fEtaH;
120   delete fEneH;
121   delete fPhiH;
122   delete fInJetH;
123   delete fFragH;
124   delete fFragLnH;
125   delete fFragrH;
126   delete fFragLnrH;
127   delete fShapeH;
128   delete fShaperH;
129 }
130
131 ////////////////////////////////////////////////////////////////////////
132
133 void AliJetControlPlots::FillHistos(AliJet *j)
134 {
135   // Fills the histograms
136   Int_t nj = j->GetNJets();
137   fNJetsH->Fill(nj,1);
138   if (nj == 0) return;
139   
140   // kinematics, occupancy and multiplicities
141   TArrayI mj = j->GetMultiplicities();
142   Int_t mjTot=0;
143   for (Int_t i=0;i<nj;i++) {
144     mjTot+=mj[i];
145     fMultH->Fill(mj[i],1);    
146     fPtH->Fill(j->GetPt(i),1);
147     fEtaH->Fill(j->GetEta(i),1);
148     fEneH->Fill(j->GetE(i),1);
149     fPhiH->Fill(j->GetPhi(i),1);
150   }
151   fInJetH->Fill(((Double_t) mjTot)/((Double_t) j->GetNinput()),1);
152   
153   // fragmentation of leading jet 
154   TArrayI inJet = j->GetInJet();
155   TArrayF etain = j->GetEtaIn();
156   TArrayF ptin = j->GetPtIn();
157   for(Int_t i=0;i<(inJet.GetSize());i++) {
158     Float_t t1 = TMath::Tan(2.0*TMath::ATan(TMath::Exp(-etain[i])));
159     Float_t ene = ptin[i]*TMath::Sqrt(1.+1./(t1*t1));
160     if (inJet[i] == 1) {
161       fFragH->Fill((Float_t) ene/(j->GetE(0)),1);
162       fFragLnH->Fill((Float_t) TMath::Log((j->GetE(0))/ene),1);
163     }
164     if (inJet[i] == -1) {
165       fFragrH->Fill((Float_t) ene/(j->GetE(0)),1);
166       fFragLnrH->Fill((Float_t) TMath::Log((j->GetE(0))/ene),1);
167     }
168   }
169   
170   // shape of leading jet 
171   // * CAREFUL: depends on number of bins and bin sizes
172   //   of shape histo. HARDWIRED at the moment *
173   //   Plot against dr NOT dr/R_jet!!!
174   TArrayF phiin = j->GetPhiIn();
175   Float_t rptS[20], rptR[20];
176   for(Int_t i=0;i<20;i++) rptS[i]=rptR[i]=0.0;
177   
178   for(Int_t i=0;i<inJet.GetSize();i++) {
179     if (inJet[i] == 1 || inJet[i] == -1) {
180       Float_t deta = etain[i] - j->GetEta(0);
181       Float_t dphi = phiin[i] - j->GetPhi(0);
182       if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
183       if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
184       Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
185       if (dr>1) continue;
186       Int_t bin = (Int_t) TMath::Floor(dr/0.05);
187       if (inJet[i] == 1) rptS[bin]+=ptin[i]/(j->GetPt(0));
188       if (inJet[i] == -1) rptR[bin]+=ptin[i]/(j->GetPt(0));       
189     }
190   }
191   
192   Float_t ptS,ptR,r;
193   ptS=ptR=0.0;
194   for (Int_t i=0;i<20;i++) {
195     ptS+=rptS[i];
196     ptR+=rptR[i];
197     r=(i+1)*0.05-0.025;
198     fShapeH->Fill(r,ptS);
199     fShaperH->Fill(r,ptR);      
200   }
201   fNJetT++;
202 }
203
204 ////////////////////////////////////////////////////////////////////////
205
206 void AliJetControlPlots::Normalize()
207 {
208 // *CAREFUL: depends on histogram number of bins 
209   if (fNJetT == 0) return;
210   fFragH->Scale(20.0/((Double_t) fNJetT));
211   fFragLnH->Scale(2.0/((Double_t) fNJetT));
212   fFragrH->Scale(20.0/((Double_t) fNJetT));
213   fFragLnrH->Scale(2.0/((Double_t) fNJetT));
214   fShapeH->Scale(1.0/((Double_t) fNJetT));
215   fShaperH->Scale(1.0/((Double_t) fNJetT));
216 }
217
218 ////////////////////////////////////////////////////////////////////////
219
220 void AliJetControlPlots::PlotHistos()
221 {
222   // style
223   gStyle->SetOptStat(kFALSE);
224   gStyle->SetOptTitle(kFALSE);
225   
226   TCanvas *c = new TCanvas("c","AliJetControlPlots",50,200,1200,700);
227   c->Divide(4,3);
228   c->cd(1);  fNJetsH->Draw("e1");
229   c->cd(2);  fMultH->Draw("e1");
230   c->cd(3);  fInJetH->Draw("e1");
231   c->cd(4);  fPtH->Draw("e1");
232   c->cd(5);  fEtaH->Draw("e1");
233   c->cd(6);  fPhiH->Draw("e1");
234   c->cd(7);  fEneH->Draw("e1");
235   c->cd(8);  fFragLnH->Draw("e1");
236   c->cd(9);  fFragLnrH->Draw("e1");
237   c->cd(10);  fShapeH->Draw("e1");
238   c->cd(11);  fShaperH->Draw("e1");
239 }
240
241 ////////////////////////////////////////////////////////////////////////
242
243 void AliJetControlPlots::SetProperties(TH1* h,const char* x, const char* y) const
244 {
245 // Sets the histogram style properties
246   h->SetMarkerStyle(20);
247   h->SetMarkerSize(.5);
248   h->SetMarkerColor(2);
249   h->SetXTitle(x);
250   h->SetYTitle(y);
251   // h->Sumw2();
252 }
253
254