1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 //---------------------------------------------------------------------
18 // FastJet v2.3.4 finder algorithm interface
20 // Author: swensy.jangal@ires.in2p3.fr
22 // Last modification: Neutral cell energy included in the jet reconstruction
24 // Author: Magali.estienne@subatech.in2p3.fr
25 //---------------------------------------------------------------------
28 #include <Riostream.h>
33 #include <TLorentzVector.h>
35 #include <TClonesArray.h>
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"
45 #include "fastjet/AreaDefinition.hh"
46 #include "fastjet/ClusterSequenceArea.hh"
47 #include "fastjet/JetDefinition.hh"
48 #include "fastjet/PseudoJet.hh"
50 // get info on how fastjet was configured
51 #include "fastjet/config.h"
53 #ifdef ENABLE_PLUGIN_SISCONE
54 #include "fastjet/SISConePlugin.hh"
57 #include<sstream> // needed for internal io
64 ClassImp(AliSISConeJetFinder)
67 //____________________________________________________________________________
69 AliSISConeJetFinder::AliSISConeJetFinder():
75 //____________________________________________________________________________
77 AliSISConeJetFinder::~AliSISConeJetFinder()
82 //______________________________________________________________________________
83 void AliSISConeJetFinder::FindJets()
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();
92 //******************************** SISCONE PLUGIN CONFIGURATION
93 // Here we look for SISCone parameters in the header and we define our plugin.
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?
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
104 fastjet::JetDefinition::Plugin * plugin;
105 plugin = new fastjet::SISConePlugin(coneRadius, overlapThreshold, nPassMax, ptProtoJetMin, caching);
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.
113 TClonesArray *lvArray = fReader->GetMomentumArray();
116 // We check if lvArray is ok
119 cout << "Could not get the momentum array" << endl;
124 Int_t nIn = lvArray->GetEntries();
126 if(nIn == 0)// nIn = Number of particles in the event
128 if (debug) cout << "entries = 0 ; Event empty !!!" << endl ;
135 // Load input vectors
136 for(Int_t i = 0; i < nIn; i++)
138 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
144 fastjet::PseudoJet inputPart(px,py,pz,en);
145 inputPart.set_user_index(i);
146 inputParticles.push_back(inputPart);
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;
159 for(Int_t i=0; i<nIn; i++)
161 AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(i);
163 if(uArray->GetUnitEnergy()>0.){
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;
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
181 } // End loop on UnitArray
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
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.
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;
203 fastjet::AreaDefinition areaDef;
205 if (areaTypeNumber < 5)
207 fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea, gridScatter, ktScatter, meanGhostKt);
208 areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
211 if (areaTypeNumber == 5)
213 Double_t effectiveRFact = header->GetEffectiveRFact();
214 fastjet::VoronoiAreaSpec ghostSpec(effectiveRFact);
215 areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
218 //********************************
219 //********************************
221 Bool_t bgMode = header->GetBGMode();// Here one choose to subtract BG or not
223 //********************************
224 //********************************
226 //********************************
227 //******************************** BG SUBTRACTION
228 if (bgMode == 1)// BG subtraction
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);
237 //***************************** BACKGROUND SUBTRACTION
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);
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;
257 // Vector of corrected jets
258 vector<fastjet::PseudoJet> corrJets = clustSeq.subtracted_jets(rho, ptMin);
260 //***************************** JETS DISPLAY
262 for (size_t j = 0; j < jets.size(); j++)
264 // If the jet is only made of ghosts, continue.
265 if (clustSeq.is_pure_ghost(jets[j]) == 1) continue;
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;
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;
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;
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("");
286 //********************************
287 //******************************** NO BG SUBTRACTION
289 if (bgMode == 0)// No BG subtraction
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);
298 //***************************** JETS DISPLAY
300 for (size_t k = 0; k < jets.size(); k++)
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;
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("");
322 //____________________________________________________________________________
324 Float_t AliSISConeJetFinder::EtaToTheta(Float_t arg)
326 // return (180./TMath::Pi())*2.*atan(exp(-arg));
327 return 2.*atan(exp(-arg));
332 //____________________________________________________________________________
334 void AliSISConeJetFinder::InitTask(TChain *tree)
337 printf("SISCone jet finder initialization ******************");
338 fReader->CreateTasks(tree);
343 //____________________________________________________________________________
345 void AliSISConeJetFinder::WriteJHeaderToFile() const
350 //____________________________________________________________________________