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