]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliSISConeJetFinder.cxx
Gsatt replaced
[u/mrichter/AliRoot.git] / JETAN / AliSISConeJetFinder.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 //---------------------------------------------------------------------
18 // FastJet v2.3.4 finder algorithm interface
19 //
20 // Author: swensy.jangal@ires.in2p3.fr 
21 //
22 // Last modification: Neutral cell energy included in the jet reconstruction
23 //
24 // Author: Magali.estienne@subatech.in2p3.fr 
25 //---------------------------------------------------------------------
26
27
28 #include <Riostream.h>
29 #include <TArrayF.h>
30 #include <TFile.h>
31 #include <TH1F.h>
32 #include <TH2F.h>
33 #include <TLorentzVector.h>
34 #include <TRandom.h>
35 #include <TClonesArray.h>
36
37 #include "AliHeader.h"
38 #include "AliJetKineReader.h"
39 #include "AliJetReader.h"
40 #include "AliJetReaderHeader.h"
41 #include "AliJetUnitArray.h"
42 #include "AliSISConeJetFinder.h"
43 #include "AliSISConeJetHeader.h"
44
45 #include "fastjet/AreaDefinition.hh"
46 #include "fastjet/ClusterSequenceArea.hh"
47 #include "fastjet/JetDefinition.hh"
48 #include "fastjet/PseudoJet.hh"
49
50 // get info on how fastjet was configured
51 #include "fastjet/config.h"
52
53 #ifdef ENABLE_PLUGIN_SISCONE
54 #include "fastjet/SISConePlugin.hh"
55 #endif
56
57 #include<sstream>  // needed for internal io
58 #include<vector> 
59 #include <cmath> 
60
61 using namespace std;
62
63
64 ClassImp(AliSISConeJetFinder)
65
66
67 //____________________________________________________________________________
68
69 AliSISConeJetFinder::AliSISConeJetFinder():
70 AliJetFinder()
71 {
72   // Constructor
73 }
74
75 //____________________________________________________________________________
76
77 AliSISConeJetFinder::~AliSISConeJetFinder()
78 {
79   // destructor
80 }
81
82 //______________________________________________________________________________
83 void AliSISConeJetFinder::FindJets()
84 {
85
86   // Pick up siscone header
87   AliSISConeJetHeader *header = (AliSISConeJetHeader*)fHeader;
88   Int_t debug  = header->GetDebug();     // debug option
89   Int_t fOpt    = fReader->GetReaderHeader()->GetDetector();
90
91
92 //******************************** SISCONE PLUGIN CONFIGURATION
93 // Here we look for SISCone parameters in the header and we define our plugin.  
94
95   Double_t coneRadius                = header->GetConeRadius();                 // cone radius
96   Double_t overlapThreshold          = header->GetOverlapThreshold();           // overlap parameter
97   Int_t    nPassMax                  = header->GetNPassMax();                   // maximum number of passes
98   Double_t ptProtoJetMin             = header->GetPtProtojetMin();              // pT min of protojets
99   Double_t caching                   = header->GetCaching();                    // do we record found cones for this set of data?
100
101 //  if (header->GetSplitMergeScale() == 0) fastjet::SISConePlugin::SplitMergeScale splitMergeScale = fastjet::SISConePlugin::SM_pttilde; // There's only one split merge scale
102 //  Double_t splitMergeStoppingScale = header->GetSplitMergeStoppingScale(); // Additional cut on pt_tilde of protojets
103
104   fastjet::JetDefinition::Plugin * plugin;
105   plugin = new fastjet::SISConePlugin(coneRadius, overlapThreshold, nPassMax, ptProtoJetMin, caching);
106
107   vector<fastjet::PseudoJet> inputParticles; 
108 //********************************  READING OF INPUT PARTICLES
109 // Here we look for px, py pz and energy of each particle that we gather in a PseudoJet object, and we put all these PseudoJet in a vector of PseudoJets : input_particles. 
110
111     if(fOpt==0) 
112       {
113         TClonesArray *lvArray = fReader->GetMomentumArray();
114
115         
116         // We check if lvArray is ok
117         if(lvArray == 0)
118           {
119             cout << "Could not get the momentum array" << endl;
120             delete plugin;
121             return;
122           }
123         
124         Int_t nIn = lvArray->GetEntries();
125         
126         if(nIn == 0)// nIn = Number of particles in the event
127           {
128             if (debug) cout << "entries = 0 ; Event empty !!!" << endl ;
129             delete plugin;
130             return;
131           }
132         
133         Float_t px,py,pz,en;
134         
135         // Load input vectors
136         for(Int_t i = 0; i < nIn; i++)
137           { 
138             TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
139             px = lv->Px();
140             py = lv->Py();
141             pz = lv->Pz();
142             en = lv->Energy();
143             
144             fastjet::PseudoJet inputPart(px,py,pz,en); 
145             inputPart.set_user_index(i);
146             inputParticles.push_back(inputPart); 
147             
148           }
149       }
150     else {
151       TClonesArray* fUnit = fReader->GetUnitArray();
152       if(fUnit == 0) { cout << "Could not get the momentum array" << endl;          delete plugin; return; }
153       Int_t         nIn = fUnit->GetEntries();
154       if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ;           delete plugin; return; }
155         // Information extracted from fUnitArray
156       // load input vectors and calculate total energy in array
157       Float_t pt,eta,phi,theta,px,py,pz,en;
158       Int_t ipart = 0;
159       for(Int_t i=0; i<nIn; i++) 
160         {
161           AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(i);
162           
163           if(uArray->GetUnitEnergy()>0.){
164             
165             // It is not necessary anymore to cut on particle pt
166             pt    = uArray->GetUnitEnergy();
167             eta   = uArray->GetUnitEta();
168             phi   = uArray->GetUnitPhi();
169             theta = EtaToTheta(eta);
170             en    = (TMath::Abs(TMath::Sin(theta)) == 0) ? pt : pt/TMath::Abs(TMath::Sin(theta));
171             px    = TMath::Cos(phi)*pt;
172             py    = TMath::Sin(phi)*pt;
173             pz    = en*TMath::TanH(eta);
174             if(debug) cout << "pt: " << pt << ", eta: " << eta << ", phi: " << phi << ", en: " << en << ", px: " << px << ", py: " << py << ", pz: " << pz << endl;
175             
176             fastjet::PseudoJet inputPart(px,py,pz,en); // create PseudoJet object
177             inputPart.set_user_index(ipart); //label the particle into Fastjet algortihm
178             inputParticles.push_back(inputPart);  // back of the input_particles vector 
179             ipart++;
180           }
181         } // End loop on UnitArray 
182     }
183   
184 //******************************** CHOICE OF JET AREA  
185 // Here we determine jets area for subtracting background later
186 // For more informations about jet areas see : The Catchment Area of Jets M. Cacciari, G. Salam and G. Soyez
187
188   Double_t ghostEtamax        = header->GetGhostEtaMax();       // maximum eta in which a ghost can be generated
189   Double_t ghostArea          = header->GetGhostArea();         // area of a ghost
190   Int_t    activeAreaRepeats  = header->GetActiveAreaRepeats(); // do we repeat area calculation?
191   Double_t gridScatter        = header->GetGridScatter();       // fractional random fluctuations of the position of the ghosts on the y-phi grid
192   Double_t ktScatter          = header->GetKtScatter();         // fractional random fluctuations of the tranverse momentum of the ghosts on the y-phi grid
193   Double_t meanGhostKt        = header->GetMeanGhostKt();       // average transverse momentum of the ghosts.
194
195   Double_t areaTypeNumber  = header->GetAreaTypeNumber();       // the number determines jet area type 
196   fastjet::AreaType areaType = fastjet::active_area;
197   if (areaTypeNumber == 1)  areaType = fastjet::active_area;                 
198   if (areaTypeNumber == 2)  areaType = fastjet::active_area_explicit_ghosts; 
199   if (areaTypeNumber == 3)  areaType = fastjet::one_ghost_passive_area;      
200   if (areaTypeNumber == 4)  areaType = fastjet::passive_area;                
201   if (areaTypeNumber == 5)  areaType = fastjet::voronoi_area;                
202
203   fastjet::AreaDefinition areaDef;
204
205   if (areaTypeNumber < 5) 
206   {
207     fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea, gridScatter, ktScatter, meanGhostKt);
208     areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
209   }
210
211   if (areaTypeNumber == 5)
212   {
213     Double_t effectiveRFact = header->GetEffectiveRFact();
214     fastjet::VoronoiAreaSpec ghostSpec(effectiveRFact);
215     areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
216   }
217
218 //********************************
219 //********************************
220
221   Bool_t bgMode = header->GetBGMode();// Here one choose to subtract BG or not
222
223 //********************************
224 //********************************
225
226 //********************************
227 //******************************** BG SUBTRACTION
228   if (bgMode == 1)// BG subtraction
229   {
230     //******************************** JETS FINDING AND EXTRACTION
231     fastjet::ClusterSequenceArea clustSeq(inputParticles, plugin, areaDef);
232     // Here we extract inclusive jets with pt > ptmin, sorted by pt 
233     Double_t ptMin = header->GetMinJetPt(); 
234     vector<fastjet::PseudoJet> inclusiveJets = clustSeq.inclusive_jets(ptMin);
235     vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets);
236   
237     //***************************** BACKGROUND SUBTRACTION
238
239     // Set the rapidity-azimuth range within which to study background
240     Double_t rapMin = header->GetRapMin();
241     Double_t rapMax = header->GetRapMax();
242     Double_t phiMin = header->GetPhiMin();
243     Double_t phiMax = header->GetPhiMax();
244     fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
245
246     // rho calculus from Cambridge/Aachen clustering (not from SISCONE as it gives too small areas)
247     Int_t algo = header->GetBGAlgorithm();
248     fastjet::JetAlgorithm algorithm = fastjet::kt_algorithm;
249     if (algo == 0) algorithm = fastjet::kt_algorithm;
250     if (algo == 1) algorithm = fastjet::cambridge_algorithm;
251     Double_t RRho = header->GetRForRho();
252     fastjet::JetDefinition jetDefForRho(algorithm, RRho);
253     fastjet::ClusterSequenceArea csForRho(inputParticles, jetDefForRho, areaDef);
254     Double_t rho = csForRho.median_pt_per_unit_area_4vector(range);
255     cout<<"rho = "<<rho<<endl;
256
257     // Vector of corrected jets
258     vector<fastjet::PseudoJet> corrJets = clustSeq.subtracted_jets(rho, ptMin);
259   
260     //***************************** JETS DISPLAY
261
262     for (size_t j = 0; j < jets.size(); j++)
263     {
264       // If the jet is only made of ghosts, continue.
265       if (clustSeq.is_pure_ghost(jets[j]) == 1) continue;
266
267       // If the correction is > jet energy px = py = pz = e = 0
268       if (corrJets[j].px() == 0 && corrJets[j].py() == 0 && corrJets[j].pz() == 0 && corrJets[j].E() == 0) continue;
269
270       if(debug){
271         cout<<"********************************** Reconstructed jet(s) (non corrected)"<<endl;
272         cout<<"Jet number "<<j+1<<" : "<<"Rapidity : "<<jets[j].rap()<<" Phi : "<<jets[j].phi()<<" pT : "<<jets[j].perp()<<" area : "<<clustSeq.area(jets[j])<<endl;
273         cout<<"e = "<<jets[j].E()<<endl;
274         
275         cout<<"********************************** Corrected jet(s)"<<endl;
276         cout<<"Jet number "<<j+1<<" : "<<"Rapidity : "<<corrJets[j].rap()<<" Phi : "<<corrJets[j].phi()<<" pT : "<<corrJets[j].perp()<<endl;
277         cout<<"e = "<<corrJets[j].E()<<endl;
278       }
279       // Go to write AOD info
280       AliAODJet aodjet (corrJets[j].px(), corrJets[j].py(), corrJets[j].pz(), corrJets[j].E());
281       if(debug) aodjet.Print("");
282       AddJet(aodjet);
283     }
284   }
285
286 //********************************
287 //******************************** NO BG SUBTRACTION
288
289   if (bgMode == 0)// No BG subtraction
290   {
291     //******************************** JETS FINDING AND EXTRACTION
292     fastjet::ClusterSequenceArea clustSeq(inputParticles, plugin, areaDef);
293     // Here we extract inclusive jets with pt > ptmin, sorted by pt 
294     Double_t ptMin = header->GetMinJetPt(); 
295     vector<fastjet::PseudoJet> inclusiveJets = clustSeq.inclusive_jets(ptMin);
296     vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets);
297
298     //***************************** JETS DISPLAY
299
300     for (size_t k = 0; k < jets.size(); k++)
301     {
302       if(debug)
303       {
304         cout<<"********************************** Reconstructed jet(s) (non corrected)"<<endl;
305         cout<<"Jet number "<<k+1<<" : "<<"Rapidity : "<<jets[k].rap()<<" Phi : "<<jets[k].phi()<<" pT : "<<jets[k].perp()<<endl;
306         cout<<"e = "<<jets[k].E()<<endl;
307       }
308
309       // Go to write AOD info
310       Double_t area = clustSeq.area(jets[k]);
311       AliAODJet aodjet (jets[k].px(), jets[k].py(), jets[k].pz(), jets[k].E());
312       aodjet.SetEffArea(area,0);
313       if(debug) aodjet.Print("");
314       AddJet(aodjet);
315     }
316   }
317
318   delete plugin;
319
320 }
321  
322 //____________________________________________________________________________
323
324 Float_t  AliSISConeJetFinder::EtaToTheta(Float_t arg)
325 {
326   //  return (180./TMath::Pi())*2.*atan(exp(-arg));
327   return 2.*atan(exp(-arg));
328
329
330 }
331
332 //____________________________________________________________________________
333
334 void AliSISConeJetFinder::InitTask(TChain *tree)
335 {
336
337   printf("SISCone jet finder initialization ******************");
338   fReader->CreateTasks(tree);
339
340 }
341
342
343 //____________________________________________________________________________
344
345 void AliSISConeJetFinder::WriteJHeaderToFile() const
346 {
347   fHeader->Write();
348 }
349
350 //____________________________________________________________________________
351