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