]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliJetFinder.cxx
Updates provided by Magali Estienne
[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 //          magali.estienne@subatech.in2p3.fr
24 //---------------------------------------------------------------------
25
26 #include <Riostream.h>
27 #include <TFile.h>
28
29 #include "AliJetFinder.h"
30 #include "AliJet.h"
31 #include "AliAODJet.h"
32 #include "AliJetControlPlots.h"
33 #include "AliLeading.h"
34 #include "AliAODEvent.h"
35 #include "AliJetUnitArray.h"
36
37 class TProcessID;
38 class TClonesArray;
39
40 ClassImp(AliJetFinder)
41
42 AliJetFinder::AliJetFinder():
43     fTreeJ(0),
44     fPlotMode(kFALSE),
45     fJets(0),
46     fGenJets(0),
47     fLeading(0),
48     fReader(0x0),
49     fHeader(0x0),
50     fAODjets(0x0),
51     fNAODjets(0),
52     fPlots(0x0),
53     fOut(0x0)
54     
55 {
56   //
57   // Constructor
58   //
59
60   fJets    = new AliJet();
61   fGenJets = new AliJet();
62   fLeading = new AliLeading();
63   fAODjets = 0;
64 }
65
66 ////////////////////////////////////////////////////////////////////////
67 AliJetFinder::~AliJetFinder()
68 {
69   //
70   // Destructor
71   //
72
73   // Reset and delete jets
74   fJets->ClearJets();
75   delete fJets;
76   fGenJets->ClearJets();
77   delete fGenJets;
78   // close file
79   if (fOut) {
80     fOut->Close();
81     fOut->Delete();
82   }
83   delete fOut;
84   // reset and delete control plots
85   if (fPlots) delete fPlots;
86 }
87
88 ////////////////////////////////////////////////////////////////////////
89 void AliJetFinder::SetOutputFile(const char */*name*/)
90 {
91   //  opens output file 
92   //  fOut = new TFile(name,"recreate");
93 }
94
95 ////////////////////////////////////////////////////////////////////////
96 void AliJetFinder::PrintJets()
97 {
98   // Print jet information
99   cout << " Jets found with jet algorithm:" << endl;
100   fJets->PrintJets();
101   cout << " Jets found by pythia:" << endl;
102   fGenJets->PrintJets();
103 }
104
105 ////////////////////////////////////////////////////////////////////////
106 void AliJetFinder::SetPlotMode(Bool_t b)
107 {
108   // Sets the plotting mode
109   fPlotMode=b;
110   if (b && !fPlots) fPlots = new AliJetControlPlots(); 
111 }
112
113 ////////////////////////////////////////////////////////////////////////
114 TTree* AliJetFinder::MakeTreeJ(char* name)
115 {
116     // Create the tree for reconstructed jets
117     fTreeJ = new TTree(name, "AliJet");
118     fTreeJ->Branch("FoundJet",   &fJets,   1000);
119     fTreeJ->Branch("GenJet",     &fGenJets,1000);
120     fTreeJ->Branch("LeadingPart",&fLeading,1000);
121     return fTreeJ;
122 }
123
124 ////////////////////////////////////////////////////////////////////////
125 void AliJetFinder::WriteRHeaderToFile()
126 {
127   // write reader header
128     AliJetReaderHeader *rh = fReader->GetReaderHeader();
129     rh->Write();
130 }
131
132 ////////////////////////////////////////////////////////////////////////
133 void AliJetFinder::Run()
134 {
135   // Do some initialization
136   Init();
137   // connect files
138   fReader->OpenInputFiles();
139
140   // write headers
141   WriteHeaders();
142   // loop over events
143   Int_t nFirst, nLast, option, debug, arrayInitialised; 
144   nFirst = fReader->GetReaderHeader()->GetFirstEvent();
145   nLast  = fReader->GetReaderHeader()->GetLastEvent();
146
147   option = fReader->GetReaderHeader()->GetDetector();
148   debug  = fReader->GetReaderHeader()->GetDebug();
149   arrayInitialised = fReader->GetArrayInitialised();
150
151   // loop over events
152   for (Int_t i=nFirst;i<nLast;i++) {
153       fReader->FillMomentumArray();
154       fLeading->FindLeading(fReader);
155       fReader->GetGenJets(fGenJets);
156
157       if (option == 0) { // TPC with fMomentumArray
158           if(debug > 1) 
159               printf("In FindJetsC() routine: find jets with fMomentumArray !!!\n");
160           FindJetsC();
161       } else {
162         if(debug > 1) printf("In FindJets() routine: find jets with fUnitArray !!!\n");
163         FindJets();
164       }
165       if (fOut) {
166           fOut->cd();
167       }
168       
169       if (fPlots) fPlots->FillHistos(fJets);
170       fLeading->Reset();
171       fGenJets->ClearJets();
172       Reset();
173   } 
174
175   // write out
176   if (fPlots) {
177       fPlots->Normalize();
178       fPlots->PlotHistos();
179   }
180   if (fOut) {
181       fOut->cd();
182       fPlots->Write();
183       fOut->Close();
184   }
185 }
186
187
188 //
189 // The following methods have been added to allow for event steering from the outside
190 //
191
192 ////////////////////////////////////////////////////////////////////////
193 void AliJetFinder::ConnectTree(TTree* tree, TObject* data)
194 {
195     // Connect the input file
196     fReader->ConnectTree(tree, data);
197 }
198
199 ////////////////////////////////////////////////////////////////////////
200 void AliJetFinder::WriteHeaders()
201 {
202     // Write the Headers
203     TFile* f = new TFile("jets_local.root", "recreate");
204     WriteRHeaderToFile();
205     WriteJHeaderToFile();
206     f->Close();
207 }
208
209 ////////////////////////////////////////////////////////////////////////
210 Bool_t AliJetFinder::ProcessEvent()
211 {
212   //
213   // Process one event
214   // Charged only jets
215   //
216
217   Bool_t ok = fReader->FillMomentumArray();
218   if (!ok) return kFALSE;
219   
220   // Leading particles
221   fLeading->FindLeading(fReader);
222   // Jets
223   FindJets(); // V1
224   //  FindJetsC(); // V2
225   
226   if (fPlots) fPlots->FillHistos(fJets);
227   fLeading->Reset();
228   fGenJets->ClearJets();
229   Reset();  
230   return kTRUE;
231 }
232
233 ////////////////////////////////////////////////////////////////////////
234 Bool_t AliJetFinder::ProcessEvent2()
235 {
236   //
237   // Process one event
238   // Charged only or charged+neutral jets
239   //
240
241   TRefArray* ref = new TRefArray();
242   Bool_t procid = kFALSE;
243   Bool_t ok = fReader->ExecTasks(procid,ref);
244
245   // Delete reference pointer  
246   if (!ok) {delete ref; return kFALSE;}
247   
248   // Leading particles
249   fLeading->FindLeading(fReader);
250   // Jets
251   FindJets();
252   
253   Int_t         nEntRef    = ref->GetEntries();
254
255   for(Int_t i=0; i<nEntRef; i++)
256     { 
257       // Reset the UnitArray content which were referenced
258       ((AliJetUnitArray*)ref->At(i))->SetUnitTrackID(0);
259       ((AliJetUnitArray*)ref->At(i))->SetUnitEnergy(0.);
260       ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag(kPtSmaller);
261       ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag2(kPtSmaller);
262       ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlag(kBad);
263       ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlagC(kTRUE,kBad);
264       ((AliJetUnitArray*)ref->At(i))->SetUnitDetectorFlag(kTpc);
265       ((AliJetUnitArray*)ref->At(i))->SetUnitFlag(kOutJet);
266       ((AliJetUnitArray*)ref->At(i))->ClearUnitTrackRef();
267
268       // Reset process ID
269       AliJetUnitArray* uA = (AliJetUnitArray*)ref->At(i);
270       uA->ResetBit(kIsReferenced);
271       uA->SetUniqueID(0);     
272     }
273
274   // Delete the reference pointer
275   ref->Delete();
276   delete ref;
277
278   if (fPlots) fPlots->FillHistos(fJets);
279   fLeading->Reset();
280   fGenJets->ClearJets();
281   Reset();
282
283   return kTRUE;
284 }
285
286 ////////////////////////////////////////////////////////////////////////
287 void AliJetFinder::FinishRun()
288 {
289     // Finish a run
290     if (fPlots) {
291         fPlots->Normalize();
292         fPlots->PlotHistos();
293     }
294     
295     if (fOut) {
296          fOut->cd();
297          if (fPlots) {
298              fPlots->Write();
299          }
300          fOut->Close();
301     }
302 }
303
304 ////////////////////////////////////////////////////////////////////////
305 void AliJetFinder::AddJet(AliAODJet p)
306 {
307 // Add new jet to the list
308   new ((*fAODjets)[fNAODjets++]) AliAODJet(p);
309 }
310
311 ////////////////////////////////////////////////////////////////////////
312 void AliJetFinder::ConnectAOD(AliAODEvent* aod)
313 {
314 // Connect to the AOD
315     fAODjets = aod->GetJets();
316 }
317
318 ////////////////////////////////////////////////////////////////////////
319 void AliJetFinder::ConnectAODNonStd(AliAODEvent* aod,const char *bname)
320 {
321
322   fAODjets = dynamic_cast<TClonesArray*>(aod->FindListObject(bname));
323   // how is this is reset? Cleared?
324 }
325