]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetFinder.cxx
Standard connection to ESD in chain.
[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     fOut = new TFile("jets.root","recreate");
107     fOut->cd();
108     fTreeJ = new TTree(name, "AliJet");
109     fTreeJ->Branch("FoundJet",   &fJets,   1000);
110     fTreeJ->Branch("GenJet",     &fGenJets,1000);
111     fTreeJ->Branch("LeadingPart",&fLeading,1000);
112     return fTreeJ;
113 }
114
115 ////////////////////////////////////////////////////////////////////////
116
117 void AliJetFinder::WriteRHeaderToFile()
118 {
119   // write reader header
120   fOut->cd();
121   AliJetReaderHeader *rh = fReader->GetReaderHeader();
122   rh->Write();
123 }
124
125 ////////////////////////////////////////////////////////////////////////
126
127
128 void AliJetFinder::Run()
129 {
130   // Do some initialization
131   Init();
132   // connect files
133   fReader->OpenInputFiles();
134
135   // write headers
136   WriteHeaders();
137   // loop over events
138   Int_t nFirst, nLast, option, debug, arrayInitialised; 
139   nFirst = fReader->GetReaderHeader()->GetFirstEvent();
140   nLast  = fReader->GetReaderHeader()->GetLastEvent();
141
142   option = fReader->GetReaderHeader()->GetDetector();
143   debug  = fReader->GetReaderHeader()->GetDebug();
144   arrayInitialised = fReader->GetArrayInitialised();
145
146   // loop over events
147   for (Int_t i=nFirst;i<nLast;i++) {
148       fReader->FillMomentumArray(i);
149       fLeading->FindLeading(fReader);
150       fReader->GetGenJets(fGenJets);
151
152       if (option == 0) { // TPC with fMomentumArray
153           if(debug > 1) 
154               printf("In FindJetsTPC() routine: find jets with fMomentumArray !!!\n");
155           FindJetsTPC();
156       } else {
157            if(debug > 1) printf("In FindJets() routine: find jets with fUnitArray !!!\n");
158            FindJets();
159       }
160       if (fOut) {
161           fOut->cd();
162           fTreeJ->Fill();
163       }
164       
165       if (fPlots) fPlots->FillHistos(fJets);
166       fLeading->Reset();
167       fGenJets->ClearJets();
168       Reset();
169   } 
170
171   // write out
172   if (fPlots) {
173       fPlots->Normalize();
174       fPlots->PlotHistos();
175   }
176   if (fOut) {
177       fOut->cd();
178       fPlots->Write();
179       fOut->Close();
180   }
181 }
182
183
184 //
185 // The following methods have been added to allow for event steering from the outside
186 //
187
188 void AliJetFinder::ConnectTree(TTree* tree, TObject* data)
189 {
190     // Connect the input file
191     fReader->ConnectTree(tree, data);
192 }
193
194 void AliJetFinder::WriteHeaders()
195 {
196     // Write the Headers
197     if (fOut) {
198         fOut->cd();
199         WriteRHeaderToFile();
200         WriteJHeaderToFile();
201     }
202 }
203
204
205 Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
206 {
207 //
208 // Process one event
209 //
210     Int_t debug  = fReader->GetReaderHeader()->GetDebug();
211     if (debug > 0) printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
212     Bool_t ok = fReader->FillMomentumArray(entry);
213     if (!ok) return kFALSE;
214     fLeading->FindLeading(fReader);
215     FindJets();
216     if (fOut) {
217         fOut->cd();
218         fTreeJ->Fill();
219     }
220     
221     if (fPlots) fPlots->FillHistos(fJets);
222     fLeading->Reset();
223     fGenJets->ClearJets();
224     Reset();  
225     return kTRUE;
226 }
227
228 void AliJetFinder::FinishRun()
229 {
230     // Finish a run
231     if (fPlots) {
232         fPlots->Normalize();
233         fPlots->PlotHistos();
234     }
235     
236     if (fOut) {
237          fOut->cd();
238          fTreeJ->Write();
239          if (fPlots) {
240              fPlots->Write();
241          }
242          fOut->Close();
243     }
244 }
245