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