]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetDistributions.cxx
New version of the jet finder (Rafael)
[u/mrichter/AliRoot.git] / JETAN / AliJetDistributions.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 // JetDistributions class 
18 // Get different basic distributions
19 // Authors: mercedes.lopez.noriega@cern.ch
20 //---------------------------------------------------------------------
21  
22 #include "AliJetDistributions.h"
23 ClassImp(AliJetDistributions)
24  
25 ////////////////////////////////////////////////////////////////////////
26 // includes
27 #include <Riostream.h>
28 #include <TH1F.h>
29 #include <TH2F.h>
30 #include <TFile.h>
31 #include <TTree.h>
32 #include <TMath.h>
33   
34 #include "AliJetProductionDataPDC2004.h"
35 #include "AliJet.h"
36 #include "AliJetKineReaderHeader.h"
37 #include "AliJetESDReaderHeader.h"
38 #include "AliLeading.h"
39   
40 ////////////////////////////////////////////////////////////////////////
41 // constructor/destructor
42   
43 AliJetDistributions::AliJetDistributions()
44 {
45   // Constructor
46   fDirectory    = 0x0;   
47   fFile         = "jets.root";   
48   fEventMin     =  0;
49   fEventMax     = -1;
50   fRunMin       =  0;
51   fRunMax       = 11;
52   fPercentage   = -1.0;
53   SetReaderHeader();
54   
55   fPart = 0;
56   fGenJ = 0;
57   fRecJ = 0;
58
59   fRetaH = 0;
60   fRphiH = 0;
61   fRptH = 0;
62   fRetaphiH = 0;
63
64   fMultH = 0;
65
66   // for Analize()
67   fPythia = kFALSE;
68   fDoPart = kTRUE;
69   fDoGenJ = kTRUE;
70   fDoRecJ = kTRUE;
71
72 }
73
74 ////////////////////////////////////////////////////////////////////////
75 // define histogrames 
76
77 void AliJetDistributions::DefineHistograms()
78 {
79   // Define histograms to be used
80   fRetaH = new TH1F("fRetaH","Reconstructed eta",140,-1.2,1.2);
81   SetProperties(fRetaH,"#eta","entries");
82   fRphiH = new TH1F("fRphiH","Reconstructed phi",18,0,2.0*TMath::Pi());
83   SetProperties(fRphiH,"#phi","entries");
84   fRptH = new TH1F("fRptH","Reconstructed pt",150,0,150);
85   SetProperties(fRptH,"p_{T} (GeV/c)","entries");
86
87   fRetaphiH = new TH2F("fRetaphiH","Reconstructed eta vs. phi",140,-1.2,1.2,18,0.,2.0*TMath::Pi());
88   SetProperties(fRetaphiH,"#eta","#phi");
89
90   fMultH = new TH1F("fMultH","Reconstructed Multiplicity",1000,0,10000);
91   SetProperties(fMultH,"Multiplicity","entries");
92 }
93
94 void AliJetDistributions::SetProperties(TH1* h,const char* x, const char* y) const
95 {
96   // Properties of histograms (x title, y title and error propagation)
97   h->SetXTitle(x);
98   h->SetYTitle(y);
99   h->Sumw2();
100 }
101  
102 ////////////////////////////////////////////////////////////////////////
103 // fill histogrames 
104
105 void AliJetDistributions::FillHistograms()
106 {
107   // Run data 
108   AliJetProductionDataPDC2004* runData = new AliJetProductionDataPDC2004();
109   
110   // Loop over runs
111   TFile* jFile = 0x0;
112   for (Int_t iRun = fRunMin; iRun <= fRunMax; iRun++) {
113     char fn[20];
114     sprintf(fn,"%s/%s.root",fDirectory,(runData->GetRunTitle(iRun)).Data());
115     jFile = new TFile(fn);
116     printf("  Analyzing run: %d %s\n", iRun,fn);        
117     
118     // Get reader header and events to be looped over
119     AliJetReaderHeader *jReaderH = (AliJetReaderHeader*)(jFile->Get(fReaderHeader));
120     if (fEventMin == -1) fEventMin =  jReaderH->GetFirstEvent();
121     if (fEventMax == -1) {
122       fEventMax =  jReaderH->GetLastEvent();
123     } else {
124       fEventMax = TMath::Min(fEventMax, jReaderH->GetLastEvent());
125     }
126
127     //AliUA1JetHeader *jH = (AliUA1JetHeader *) (jFile->Get("AliUA1JetHeader"));
128     //jH->PrintParameters();
129
130     
131     // Loop over events
132     for (Int_t i = fEventMin; i < fEventMax; i++) {
133       //printf("  Analyzing run: %d  Event %d / %d \n",iRun, i, fEventMax);
134
135       // Get next tree with AliJet
136       char nameT[100];
137       sprintf(nameT, "TreeJ%d",i);
138       TTree *jetT =(TTree *)(jFile->Get(nameT));
139       if (fDoRecJ) jetT->SetBranchAddress("FoundJet",    &fRecJ);
140       if (fDoGenJ) jetT->SetBranchAddress("GenJet",      &fGenJ);
141       if (fDoPart) jetT->SetBranchAddress("LeadingPart", &fPart);
142       
143       jetT->GetEntry(0);
144
145       FillDistributions(fRecJ);
146
147       delete jetT;
148     } // end loop over events in one file
149     if (jFile) jFile->Close();
150     delete jFile;
151   } // end loop over files
152 }
153
154 void AliJetDistributions::FillDistributions(AliJet *j)
155 {
156   // Fill histrograms
157   TArrayI inJet = j->GetInJet();
158   TArrayF etain = j->GetEtaIn();
159   TArrayF ptin = j->GetPtIn();
160   TArrayF phiin = j->GetPhiIn();
161
162   fMultH->Fill(inJet.GetSize(),1);
163   for(Int_t i=0;i<inJet.GetSize();i++) {
164     fRetaH->Fill(etain[i],1);
165     fRphiH->Fill(phiin[i],1);
166     fRptH->Fill(ptin[i],1);
167     fRetaphiH->Fill(etain[i],phiin[i],1);
168   }
169 }
170
171 ////////////////////////////////////////////////////////////////////////
172 // Plot histogrames 
173
174 void AliJetDistributions::PlotHistograms()
175 {
176   // improved!
177   fMultH->Draw();
178   fRetaH->Draw();
179   fRphiH->Draw();
180   fRetaphiH->Draw();
181   fRptH->Draw();
182
183 }
184
185
186 ////////////////////////////////////////////////////////////////////////
187 // Save histogrames 
188
189 void AliJetDistributions::SaveHistograms()
190 {
191   // Save histogramas
192   TFile *fOut = new TFile(fFile,"recreate");
193   fOut->cd();
194   fMultH->Write();
195   fRetaphiH->Write();
196   fRetaH->Write();
197   fRphiH->Write();
198   fRptH->Write();
199   fOut->Close();
200 }
201
202 // main Analysis function
203
204 void AliJetDistributions::Analyze()
205 {
206   // Do the analysis
207   DefineHistograms();
208   FillHistograms();
209   PlotHistograms();
210   SaveHistograms();
211 }
212
213