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