Fork me on GitHub

source: git/modules/FastJetFinder.cc@ 7c0fcd5

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 7c0fcd5 was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

delete duplicate license file and prepend GPLv3 header to all source code files

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