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