]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetFinder.cxx
67d76bec640512ea3ef851d1e3baae34c6c08a19
[u/mrichter/AliRoot.git] / JETAN / AliJetFinder.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 /* $Id$ */
17
18 //---------------------------------------------------------------------
19 // Jet finder base class
20 // manages the search for jets 
21 // Authors: jgcn@mda.cinvestav.mx
22 //          andreas.morsch@cern.ch
23 //---------------------------------------------------------------------
24
25 #include <Riostream.h>
26 #include <TFile.h>
27 #include "AliJetFinder.h"
28 #include "AliJet.h"
29 #include "AliJetReader.h"
30 #include "AliJetReaderHeader.h"
31 #include "AliJetControlPlots.h"
32 #include "AliLeading.h"
33
34 ClassImp(AliJetFinder)
35
36 AliJetFinder::AliJetFinder():
37     fTreeJ(0),
38     fPlotMode(kFALSE),
39     fJets(0),
40     fGenJets(0),
41     fLeading(0),
42     fReader(0x0),
43     fHeader(0x0),
44     fPlots(0x0),
45     fOut(0x0)
46     
47 {
48   // Constructor
49   fJets    = new AliJet();
50   fGenJets = new AliJet();
51   fLeading = new AliLeading();
52 }
53
54 ////////////////////////////////////////////////////////////////////////
55
56 AliJetFinder::~AliJetFinder()
57 {
58   // destructor
59   // here reset and delete jets
60   fJets->ClearJets();
61   delete fJets;
62   fGenJets->ClearJets();
63   delete fGenJets;
64   // close file
65   if (fOut) {
66     fOut->Close();
67     fOut->Delete();
68   }
69   delete fOut;
70   // reset and delete control plots
71   if (fPlots) delete fPlots;
72 }
73
74 ////////////////////////////////////////////////////////////////////////
75
76 void AliJetFinder::SetOutputFile(const char */*name*/)
77 {
78   //  opens output file 
79   //  fOut = new TFile(name,"recreate");
80 }
81
82 ////////////////////////////////////////////////////////////////////////
83
84 void AliJetFinder::PrintJets()
85 {
86   // Print jet information
87   cout << " Jets found with jet algorithm:" << endl;
88   fJets->PrintJets();
89   cout << " Jets found by pythia:" << endl;
90   fGenJets->PrintJets();
91 }
92
93 ////////////////////////////////////////////////////////////////////////
94
95 void AliJetFinder::SetPlotMode(Bool_t b)
96 {
97   // Sets the plotting mode
98   fPlotMode=b;
99   if (b && !fPlots) fPlots = new AliJetControlPlots(); 
100 }
101
102 ////////////////////////////////////////////////////////////////////////
103 TTree* AliJetFinder::MakeTreeJ(char* name)
104 {
105     // Create the tree for reconstructed jets
106     fTreeJ = new TTree(name, "AliJet");
107     fTreeJ->Branch("FoundJet",   &fJets,   1000);
108     fTreeJ->Branch("GenJet",     &fGenJets,1000);
109     fTreeJ->Branch("LeadingPart",&fLeading,1000);
110     return fTreeJ;
111 }
112
113 ////////////////////////////////////////////////////////////////////////
114
115 void AliJetFinder::WriteRHeaderToFile()
116 {
117   // write reader header
118     AliJetReaderHeader *rh = fReader->GetReaderHeader();
119     rh->Write();
120 }
121
122 ////////////////////////////////////////////////////////////////////////
123
124
125 void AliJetFinder::Run()
126 {
127   // Do some initialization
128   Init();
129   // connect files
130   fReader->OpenInputFiles();
131
132   // write headers
133   WriteHeaders();
134   // loop over events
135   Int_t nFirst, nLast, option, debug, arrayInitialised; 
136   nFirst = fReader->GetReaderHeader()->GetFirstEvent();
137   nLast  = fReader->GetReaderHeader()->GetLastEvent();
138
139   option = fReader->GetReaderHeader()->GetDetector();
140   debug  = fReader->GetReaderHeader()->GetDebug();
141   arrayInitialised = fReader->GetArrayInitialised();
142
143   // loop over events
144   for (Int_t i=nFirst;i<nLast;i++) {
145       fReader->FillMomentumArray(i);
146       fLeading->FindLeading(fReader);
147       fReader->GetGenJets(fGenJets);
148
149       if (option == 0) { // TPC with fMomentumArray
150           if(debug > 1) 
151               printf("In FindJetsTPC() routine: find jets with fMomentumArray !!!\n");
152           FindJetsTPC();
153       } else {
154            if(debug > 1) printf("In FindJets() routine: find jets with fUnitArray !!!\n");
155            FindJets();
156       }
157       if (fOut) {
158           fOut->cd();
159           fTreeJ->Fill();
160       }
161       
162       if (fPlots) fPlots->FillHistos(fJets);
163       fLeading->Reset();
164       fGenJets->ClearJets();
165       Reset();
166   } 
167
168   // write out
169   if (fPlots) {
170       fPlots->Normalize();
171       fPlots->PlotHistos();
172   }
173   if (fOut) {
174       fOut->cd();
175       fPlots->Write();
176       fOut->Close();
177   }
178 }
179
180
181 //
182 // The following methods have been added to allow for event steering from the outside
183 //
184
185 void AliJetFinder::ConnectTree(TTree* tree, TObject* data)
186 {
187     // Connect the input file
188     fReader->ConnectTree(tree, data);
189 }
190
191 void AliJetFinder::WriteHeaders()
192 {
193     // Write the Headers
194     TFile* f = new TFile("jets_local.root", "recreate");
195     WriteRHeaderToFile();
196     WriteJHeaderToFile();
197     f->Close();
198 }
199
200
201 Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
202 {
203 //
204 // Process one event
205 //
206     Int_t debug  = fReader->GetReaderHeader()->GetDebug();
207     if (debug > 0) printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
208     
209     Bool_t ok = fReader->FillMomentumArray(entry);
210     if (!ok) return kFALSE;
211
212     // Leading particles
213     fLeading->FindLeading(fReader);
214     // Jets
215     FindJets();
216     // Fill the tree
217     fTreeJ->Fill();
218
219     if (fPlots) fPlots->FillHistos(fJets);
220     fLeading->Reset();
221     fGenJets->ClearJets();
222     Reset();  
223     return kTRUE;
224 }
225
226 void AliJetFinder::FinishRun()
227 {
228     // Finish a run
229     if (fPlots) {
230         fPlots->Normalize();
231         fPlots->PlotHistos();
232     }
233     
234     if (fOut) {
235          fOut->cd();
236          fTreeJ->Write();
237          if (fPlots) {
238              fPlots->Write();
239          }
240          fOut->Close();
241     }
242 }
243