Fork me on GitHub

source: git/modules/FastJetFinder.cc@ 2b900b8

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 2b900b8 was df04eb1, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

remove unused variables from FastJetFinder

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