Fork me on GitHub

source: git/modules/FastJetFinder.cc@ ab16005

ImprovedOutputFile Timing dual_readout llp
Last change on this file since ab16005 was 9687203, checked in by mselvaggi <mselvaggi@…>, 11 years ago

added nsubjettiness

  • Property mode set to 100644
File size: 11.2 KB
Line 
1
2/** \class FastJetFinder
3 *
4 * Finds jets using FastJet library.
5 *
6 * $Date$
7 * $Revision$
8 *
9 *
10 * \author P. Demin - UCL, Louvain-la-Neuve
11 *
12 */
13
14#include "modules/FastJetFinder.h"
15
16#include "classes/DelphesClasses.h"
17#include "classes/DelphesFactory.h"
18#include "classes/DelphesFormula.h"
19
20#include "ExRootAnalysis/ExRootResult.h"
21#include "ExRootAnalysis/ExRootFilter.h"
22#include "ExRootAnalysis/ExRootClassifier.h"
23
24#include "TMath.h"
25#include "TString.h"
26#include "TFormula.h"
27#include "TRandom3.h"
28#include "TObjArray.h"
29#include "TDatabasePDG.h"
30#include "TLorentzVector.h"
31
32#include <algorithm>
33#include <stdexcept>
34#include <iostream>
35#include <sstream>
36#include <vector>
37
38#include "fastjet/PseudoJet.hh"
39#include "fastjet/JetDefinition.hh"
40#include "fastjet/ClusterSequence.hh"
41#include "fastjet/Selector.hh"
42#include "fastjet/ClusterSequenceArea.hh"
43#include "fastjet/tools/JetMedianBackgroundEstimator.hh"
44
45#include "fastjet/plugins/SISCone/fastjet/SISConePlugin.hh"
46#include "fastjet/plugins/CDFCones/fastjet/CDFMidPointPlugin.hh"
47#include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh"
48
49#include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh"
50#include "fastjet/contribs/Nsubjettiness/Njettiness.hh"
51#include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh"
52#include "fastjet/contribs/Nsubjettiness/WinnerTakeAllRecombiner.hh"
53
54using namespace std;
55using namespace fastjet;
56using namespace fastjet::contrib;
57
58
59//------------------------------------------------------------------------------
60
61FastJetFinder::FastJetFinder() :
62 fPlugin(0), fDefinition(0), fAreaDefinition(0), fItInputArray(0), fRecomb(0), fNjettinessPlugin(0)
63{
64
65}
66
67//------------------------------------------------------------------------------
68
69FastJetFinder::~FastJetFinder()
70{
71
72}
73
74//------------------------------------------------------------------------------
75
76void FastJetFinder::Init()
77{
78
79 JetDefinition::Plugin *plugin = NULL;
80 JetDefinition::Recombiner *recomb = NULL;
81 NjettinessPlugin *njet_plugin = NULL;
82
83 // read eta ranges
84
85 ExRootConfParam param = GetParam("RhoEtaRange");
86 Long_t i, size;
87
88 fEtaRangeMap.clear();
89 size = param.GetSize();
90 for(i = 0; i < size/2; ++i)
91 {
92 fEtaRangeMap[param[i*2].GetDouble()] = param[i*2 + 1].GetDouble();
93 }
94
95 // define algorithm
96
97 fJetAlgorithm = GetInt("JetAlgorithm", 6);
98 fParameterR = GetDouble("ParameterR", 0.5);
99
100 fConeRadius = GetDouble("ConeRadius", 0.5);
101 fSeedThreshold = GetDouble("SeedThreshold", 1.0);
102 fConeAreaFraction = GetDouble("ConeAreaFraction", 1.0);
103 fMaxIterations = GetInt("MaxIterations", 100);
104 fMaxPairSize = GetInt("MaxPairSize", 2);
105 fIratch = GetInt("Iratch", 1);
106 fAdjacencyCut = GetInt("AdjacencyCut", 2);
107 fOverlapThreshold = GetDouble("OverlapThreshold", 0.75);
108
109 fJetPTMin = GetDouble("JetPTMin", 10.0);
110
111 //-- N(sub)jettiness parameters --
112
113 fComputeNsubjettiness = GetBool("ComputeNsubjettiness", false);
114 fBeta = GetDouble("Beta", 1.0);
115 fAxisMode = GetInt("AxisMode", 1);
116 fRcutOff = GetDouble("RcutOff", 0.8); //used only if Njettiness is used as jet clustering algo (case 8)
117 fN = GetInt("N", 2); //used only if Njettiness is used as jet clustering algo (case 8)
118
119 // --- Jet Area Parameters ---
120 fAreaAlgorithm = GetInt("AreaAlgorithm", 0);
121 fComputeRho = GetBool("ComputeRho", false);
122
123 // - ghost based areas -
124 fGhostEtaMax = GetDouble("GhostEtaMax", 5.0);
125 fRepeat = GetInt("Repeat", 1);
126 fGhostArea = GetDouble("GhostArea", 0.01);
127 fGridScatter = GetDouble("GridScatter", 1.0);
128 fPtScatter = GetDouble("PtScatter", 0.1);
129 fMeanGhostPt = GetDouble("MeanGhostPt", 1.0E-100);
130
131 // - voronoi based areas -
132 fEffectiveRfact = GetDouble("EffectiveRfact", 1.0);
133
134 switch(fAreaAlgorithm)
135 {
136 case 1:
137 fAreaDefinition = new fastjet::AreaDefinition(active_area_explicit_ghosts, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
138 break;
139 case 2:
140 fAreaDefinition = new fastjet::AreaDefinition(one_ghost_passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
141 break;
142 case 3:
143 fAreaDefinition = new fastjet::AreaDefinition(passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
144 break;
145 case 4:
146 fAreaDefinition = new fastjet::AreaDefinition(VoronoiAreaSpec(fEffectiveRfact));
147 break;
148 case 5:
149 fAreaDefinition = new fastjet::AreaDefinition(active_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
150 break;
151 default:
152 case 0:
153 fAreaDefinition = 0;
154 break;
155 }
156
157 switch(fJetAlgorithm)
158 {
159 case 1:
160 plugin = new fastjet::CDFJetCluPlugin(fSeedThreshold, fConeRadius, fAdjacencyCut, fMaxIterations, fIratch, fOverlapThreshold);
161 fDefinition = new fastjet::JetDefinition(plugin);
162 break;
163 case 2:
164 plugin = new fastjet::CDFMidPointPlugin(fSeedThreshold, fConeRadius, fConeAreaFraction, fMaxPairSize, fMaxIterations, fOverlapThreshold);
165 fDefinition = new fastjet::JetDefinition(plugin);
166 break;
167 case 3:
168 plugin = new fastjet::SISConePlugin(fConeRadius, fOverlapThreshold, fMaxIterations, fJetPTMin);
169 fDefinition = new fastjet::JetDefinition(plugin);
170 break;
171 case 4:
172 fDefinition = new fastjet::JetDefinition(fastjet::kt_algorithm, fParameterR);
173 break;
174 case 5:
175 fDefinition = new fastjet::JetDefinition(fastjet::cambridge_algorithm, fParameterR);
176 break;
177 default:
178 case 6:
179 fDefinition = new fastjet::JetDefinition(fastjet::antikt_algorithm, fParameterR);
180 break;
181 case 7:
182 recomb = new fastjet::contrib::WinnerTakeAllRecombiner();
183 fDefinition = new fastjet::JetDefinition(fastjet::antikt_algorithm, fParameterR, recomb, Best);
184 break;
185 case 8:
186 njet_plugin = new fastjet::contrib::NjettinessPlugin(fN, Njettiness::wta_kt_axes, Njettiness::unnormalized_cutoff_measure, fBeta, fRcutOff);
187 fDefinition = new fastjet::JetDefinition(njet_plugin);
188 break;
189 }
190
191
192 fPlugin = plugin;
193 fRecomb = recomb;
194 fNjettinessPlugin = njet_plugin;
195
196 ClusterSequence::print_banner();
197
198 // import input array
199
200 fInputArray = ImportArray(GetString("InputArray", "Calorimeter/towers"));
201 fItInputArray = fInputArray->MakeIterator();
202
203 // create output arrays
204
205 fOutputArray = ExportArray(GetString("OutputArray", "jets"));
206 fRhoOutputArray = ExportArray(GetString("RhoOutputArray", "rho"));
207}
208
209//------------------------------------------------------------------------------
210
211void FastJetFinder::Finish()
212{
213 if(fItInputArray) delete fItInputArray;
214 if(fDefinition) delete fDefinition;
215 if(fAreaDefinition) delete fAreaDefinition;
216 if(fPlugin) delete static_cast<JetDefinition::Plugin*>(fPlugin);
217 if(fRecomb) delete static_cast<JetDefinition::Recombiner*>(fRecomb);
218 if(fNjettinessPlugin) delete static_cast<JetDefinition::Plugin*>(fNjettinessPlugin);
219}
220
221//------------------------------------------------------------------------------
222
223void FastJetFinder::Process()
224{
225 Candidate *candidate, *constituent;
226 TLorentzVector momentum;
227 Double_t deta, dphi, detaMax, dphiMax;
228 Double_t time, weightTime, avTime;
229 Int_t number;
230 Double_t rho = 0;
231 PseudoJet jet, area;
232 vector<PseudoJet> inputList, outputList;
233 ClusterSequence *sequence;
234 map< Double_t, Double_t >::iterator itEtaRangeMap;
235
236 DelphesFactory *factory = GetFactory();
237
238 inputList.clear();
239
240 // loop over input objects
241 fItInputArray->Reset();
242 number = 0;
243 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
244 {
245 momentum = candidate->Momentum;
246 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
247 jet.set_user_index(number);
248 inputList.push_back(jet);
249 ++number;
250 }
251
252 // construct jets
253 if(fAreaDefinition)
254 {
255 sequence = new ClusterSequenceArea(inputList, *fDefinition, *fAreaDefinition);
256 }
257 else
258 {
259 sequence = new ClusterSequence(inputList, *fDefinition);
260 }
261
262 // compute rho and store it
263 if(fComputeRho && fAreaDefinition)
264 {
265 for(itEtaRangeMap = fEtaRangeMap.begin(); itEtaRangeMap != fEtaRangeMap.end(); ++itEtaRangeMap)
266 {
267 Selector select_rapidity = SelectorAbsRapRange(itEtaRangeMap->first, itEtaRangeMap->second);
268 JetMedianBackgroundEstimator estimator(select_rapidity, *fDefinition, *fAreaDefinition);
269 estimator.set_particles(inputList);
270 rho = estimator.rho();
271
272 candidate = factory->NewCandidate();
273 candidate->Momentum.SetPtEtaPhiE(rho, 0.0, 0.0, rho);
274 candidate->Edges[0] = itEtaRangeMap->first;
275 candidate->Edges[1] = itEtaRangeMap->second;
276 fRhoOutputArray->Add(candidate);
277 }
278 }
279
280 outputList.clear();
281 outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin));
282
283
284 // loop over all jets and export them
285 detaMax = 0.0;
286 dphiMax = 0.0;
287 vector<PseudoJet>::iterator itInputList, itOutputList;
288 for(itOutputList = outputList.begin(); itOutputList != outputList.end(); ++itOutputList)
289 {
290 momentum.SetPxPyPzE(itOutputList->px(), itOutputList->py(), itOutputList->pz(), itOutputList->E());
291 area.reset(0.0, 0.0, 0.0, 0.0);
292 if(fAreaDefinition) area = itOutputList->area_4vector();
293
294 candidate = factory->NewCandidate();
295
296 time=0;
297 weightTime=0;
298
299 inputList.clear();
300 inputList = sequence->constituents(*itOutputList);
301
302 for(itInputList = inputList.begin(); itInputList != inputList.end(); ++itInputList)
303 {
304 constituent = static_cast<Candidate*>(fInputArray->At(itInputList->user_index()));
305
306 deta = TMath::Abs(momentum.Eta() - constituent->Momentum.Eta());
307 dphi = TMath::Abs(momentum.DeltaPhi(constituent->Momentum));
308 if(deta > detaMax) detaMax = deta;
309 if(dphi > dphiMax) dphiMax = dphi;
310
311 time += TMath::Sqrt(constituent->Momentum.E())*(constituent->Position.T());
312 weightTime += TMath::Sqrt(constituent->Momentum.E());
313
314 candidate->AddCandidate(constituent);
315 }
316
317 avTime = time/weightTime;
318
319 candidate->Momentum = momentum;
320 candidate->Position.SetT(avTime);
321 candidate->Area.SetPxPyPzE(area.px(), area.py(), area.pz(), area.E());
322
323 candidate->DeltaEta = detaMax;
324 candidate->DeltaPhi = dphiMax;
325
326 // --- compute N-subjettiness with N = 1,2,3,4,5 ----
327
328 if(fComputeNsubjettiness)
329 {
330 Njettiness::AxesMode axisMode;
331
332 if (fAxisMode == 1) axisMode = Njettiness::wta_kt_axes;
333 if (fAxisMode == 2) axisMode = Njettiness::onepass_wta_kt_axes;
334 if (fAxisMode == 3) axisMode = Njettiness::kt_axes;
335 if (fAxisMode == 4) axisMode = Njettiness::onepass_kt_axes;
336
337 Njettiness::MeasureMode measureMode = Njettiness::unnormalized_measure;
338
339 Nsubjettiness nSub1(1, axisMode, measureMode, fBeta);
340 Nsubjettiness nSub2(2, axisMode, measureMode, fBeta);
341 Nsubjettiness nSub3(3, axisMode, measureMode, fBeta);
342 Nsubjettiness nSub4(4, axisMode, measureMode, fBeta);
343 Nsubjettiness nSub5(5, axisMode, measureMode, fBeta);
344
345 candidate -> Tau1 = nSub1(*itOutputList);
346 candidate -> Tau2 = nSub2(*itOutputList);
347 candidate -> Tau3 = nSub3(*itOutputList);
348 candidate -> Tau4 = nSub4(*itOutputList);
349 candidate -> Tau5 = nSub5(*itOutputList);
350 }
351
352
353 fOutputArray->Add(candidate);
354 }
355 delete sequence;
356}
Note: See TracBrowser for help on using the repository browser.