]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetFinder.cxx
Protection added.
[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 //---------------------------------------------------------------------
17 // Jet finder base class
18 // manages the search for jets 
19 // Authors: jgcn@mda.cinvestav.mx
20 //          andreas.morsch@cern.ch
21 //---------------------------------------------------------------------
22
23 #include <Riostream.h>
24 #include <TFile.h>
25 #include "AliJetFinder.h"
26 #include "AliJet.h"
27 #include "AliJetReader.h"
28 #include "AliJetReaderHeader.h"
29 #include "AliJetControlPlots.h"
30 #include "AliLeading.h"
31
32 ClassImp(AliJetFinder)
33
34 AliJetFinder::AliJetFinder():
35     fPlotMode(kFALSE),
36     fJets(0),
37     fGenJets(0),
38     fLeading(0),
39     fReader(0x0),
40     fPlots(0x0),
41     fOut(0x0)
42     
43 {
44   // Constructor
45   fJets    = new AliJet();
46   fGenJets = new AliJet();
47   fLeading = new AliLeading();
48 }
49
50 ////////////////////////////////////////////////////////////////////////
51
52 AliJetFinder::~AliJetFinder()
53 {
54   // destructor
55   // here reset and delete jets
56   fJets->ClearJets();
57   delete fJets;
58   fGenJets->ClearJets();
59   delete fGenJets;
60   // close file
61   if (fOut) {
62     fOut->Close();
63     fOut->Delete();
64   }
65   delete fOut;
66   // reset and delete control plots
67   if (fPlots) delete fPlots;
68 }
69
70 ////////////////////////////////////////////////////////////////////////
71
72 void AliJetFinder::SetOutputFile(const char *name)
73 {
74   // opens output file 
75   fOut = new TFile(name,"recreate");
76 }
77
78 ////////////////////////////////////////////////////////////////////////
79
80 void AliJetFinder::PrintJets()
81 {
82   // Print jet information
83   cout << " Jets found with jet algorithm:" << endl;
84   fJets->PrintJets();
85   cout << " Jets found by pythia:" << endl;
86   fGenJets->PrintJets();
87 }
88
89 ////////////////////////////////////////////////////////////////////////
90
91 void AliJetFinder::SetPlotMode(Bool_t b)
92 {
93   // Sets the plotting mode
94   fPlotMode=b;
95   if (b && !fPlots) fPlots = new AliJetControlPlots(); 
96 }
97
98 ////////////////////////////////////////////////////////////////////////
99
100 void AliJetFinder::WriteJetsToFile(Int_t i)
101 {
102   // Writes the jets to file
103   fOut->cd();
104   char hname[30];
105   sprintf(hname,"TreeJ%d",i);
106   TTree* jetT = new TTree(hname,"AliJet");
107   jetT->Branch("FoundJet",&fJets,1000);
108   jetT->Branch("GenJet",&fGenJets,1000);
109   jetT->Branch("LeadingPart",&fLeading,1000);
110   jetT->Fill();
111   jetT->Write(hname);
112   delete jetT;
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) WriteJetsToFile(i);
161       if (fPlots) fPlots->FillHistos(fJets);
162       fLeading->Reset();
163       fGenJets->ClearJets();
164       Reset();
165   } 
166
167   // write out
168   if (fPlots) {
169       fPlots->Normalize();
170       fPlots->PlotHistos();
171   }
172   if (fOut) {
173       fOut->cd();
174       fPlots->Write();
175       fOut->Close();
176   }
177 }
178
179
180 //
181 // The following methods have been added to allow for event steering from the outside
182 //
183
184 void AliJetFinder::ConnectTree(TTree* tree)
185 {
186     // Connect the input file
187     fReader->ConnectTree(tree);
188 }
189
190 void AliJetFinder::WriteHeaders()
191 {
192     // Write the Headers
193     if (fOut) {
194         fOut->cd();
195         WriteRHeaderToFile();
196         WriteJHeaderToFile();
197     }
198 }
199
200
201 Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
202 {
203 //
204 // Process one event
205 //
206     printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
207     Bool_t ok = fReader->FillMomentumArray(entry);
208     if (!ok) return kFALSE;
209     fLeading->FindLeading(fReader);
210     FindJets();
211     if (fOut)   WriteJetsToFile(entry);
212     if (fPlots) fPlots->FillHistos(fJets);
213     fLeading->Reset();
214     fGenJets->ClearJets();
215     Reset();  
216     return kTRUE;
217 }
218
219 void AliJetFinder::FinishRun()
220 {
221     // Finish a run
222   if (fPlots) {
223       fPlots->Normalize();
224       fPlots->PlotHistos();
225       if (fOut) {
226           fOut->cd();
227           fPlots->Write();
228           fOut->Close();
229       }   
230   } else {
231       if (fOut) fOut->Close();
232   }
233 }
234