Fork me on GitHub

source: git/modules/FastJetFinder.cc@ e9a3d17

ImprovedOutputFile Timing dual_readout llp
Last change on this file since e9a3d17 was 5a69ba1d, checked in by Ulrike Schnoor <schnooru@…>, 7 years ago

Merge branch 'master' of github.com:delphes/delphes

  • Property mode set to 100644
File size: 17.5 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"
[973b92a]66#include "fastjet/contribs/Nsubjettiness/ExtraRecombiners.hh"
[9687203]67
[f319c1d]68#include "fastjet/contribs/ValenciaPlugin/ValenciaPlugin.hh"
69
[de6d698]70#include "fastjet/tools/Filter.hh"
71#include "fastjet/tools/Pruner.hh"
72#include "fastjet/contribs/RecursiveTools/SoftDrop.hh"
73
[d7d2da3]74using namespace std;
75using namespace fastjet;
[9687203]76using namespace fastjet::contrib;
77
[d7d2da3]78
79//------------------------------------------------------------------------------
80
81FastJetFinder::FastJetFinder() :
[e9c0d73]82 fPlugin(0), fRecomb(0), fAxesDef(0), fMeasureDef(0), fNjettinessPlugin(0), fValenciaPlugin(0),
83 fDefinition(0), fAreaDefinition(0), fItInputArray(0)
[d7d2da3]84{
85
86}
87
88//------------------------------------------------------------------------------
89
90FastJetFinder::~FastJetFinder()
91{
92
93}
94
95//------------------------------------------------------------------------------
96
97void FastJetFinder::Init()
98{
[c6321ad]99 JetDefinition::Plugin *plugin = 0;
100 JetDefinition::Recombiner *recomb = 0;
[7278220]101 ExRootConfParam param;
[8336b6e]102 Long_t i, size;
[7278220]103 Double_t etaMin, etaMax;
104 TEstimatorStruct estimatorStruct;
[8336b6e]105
[d7d2da3]106 // define algorithm
107
108 fJetAlgorithm = GetInt("JetAlgorithm", 6);
109 fParameterR = GetDouble("ParameterR", 0.5);
110
111 fConeRadius = GetDouble("ConeRadius", 0.5);
112 fSeedThreshold = GetDouble("SeedThreshold", 1.0);
113 fConeAreaFraction = GetDouble("ConeAreaFraction", 1.0);
114 fMaxIterations = GetInt("MaxIterations", 100);
115 fMaxPairSize = GetInt("MaxPairSize", 2);
116 fIratch = GetInt("Iratch", 1);
[3ccc8586]117 fAdjacencyCut = GetInt("AdjacencyCut", 2);
[d7d2da3]118 fOverlapThreshold = GetDouble("OverlapThreshold", 0.75);
119
120 fJetPTMin = GetDouble("JetPTMin", 10.0);
121
[f319c1d]122
[9687203]123 //-- N(sub)jettiness parameters --
[e4c3fef]124
[9687203]125 fComputeNsubjettiness = GetBool("ComputeNsubjettiness", false);
126 fBeta = GetDouble("Beta", 1.0);
127 fAxisMode = GetInt("AxisMode", 1);
[e4c3fef]128 fRcutOff = GetDouble("RcutOff", 0.8); // used only if Njettiness is used as jet clustering algo (case 8)
129 fN = GetInt("N", 2); // used only if Njettiness is used as jet clustering algo (case 8)
[f319c1d]130
131 //-- Exclusive clustering for e+e- collisions --
132
133 fNJets = GetInt("NJets",2);
134 fExclusiveClustering = GetBool("ExclusiveClustering", false);
135
136 //-- Valencia Linear Collider algorithm
137 fGamma = GetDouble("Gamma", 1.0);
138 //fBeta parameter see above
139
[919f6a6]140 fMeasureDef = new NormalizedMeasure(fBeta, fParameterR);
[8abab33]141
142 switch(fAxisMode)
143 {
144 default:
145 case 1:
146 fAxesDef = new WTA_KT_Axes();
147 break;
148 case 2:
149 fAxesDef = new OnePass_WTA_KT_Axes();
150 break;
151 case 3:
152 fAxesDef = new KT_Axes();
153 break;
154 case 4:
155 fAxesDef = new OnePass_KT_Axes();
156 }
[e4c3fef]157
[de6d698]158 //-- Trimming parameters --
159
160 fComputeTrimming = GetBool("ComputeTrimming", false);
161 fRTrim = GetDouble("RTrim", 0.2);
162 fPtFracTrim = GetDouble("PtFracTrim", 0.05);
163
164
165 //-- Pruning parameters --
166
167 fComputePruning = GetBool("ComputePruning", false);
168 fZcutPrun = GetDouble("ZcutPrun", 0.1);
169 fRcutPrun = GetDouble("RcutPrun", 0.5);
170 fRPrun = GetDouble("RPrun", 0.8);
171
172 //-- SoftDrop parameters --
173
174 fComputeSoftDrop = GetBool("ComputeSoftDrop", false);
175 fBetaSoftDrop = GetDouble("BetaSoftDrop", 0.0);
176 fSymmetryCutSoftDrop = GetDouble("SymmetryCutSoftDrop", 0.1);
177 fR0SoftDrop= GetDouble("R0SoftDrop=", 0.8);
178
179
[d7d2da3]180 // --- Jet Area Parameters ---
181 fAreaAlgorithm = GetInt("AreaAlgorithm", 0);
182 fComputeRho = GetBool("ComputeRho", false);
[e4c3fef]183
[d7d2da3]184 // - ghost based areas -
185 fGhostEtaMax = GetDouble("GhostEtaMax", 5.0);
186 fRepeat = GetInt("Repeat", 1);
187 fGhostArea = GetDouble("GhostArea", 0.01);
188 fGridScatter = GetDouble("GridScatter", 1.0);
189 fPtScatter = GetDouble("PtScatter", 0.1);
190 fMeanGhostPt = GetDouble("MeanGhostPt", 1.0E-100);
[e4c3fef]191
[d7d2da3]192 // - voronoi based areas -
193 fEffectiveRfact = GetDouble("EffectiveRfact", 1.0);
194
195 switch(fAreaAlgorithm)
196 {
197 case 1:
[c6321ad]198 fAreaDefinition = new AreaDefinition(active_area_explicit_ghosts, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
[d7d2da3]199 break;
200 case 2:
[c6321ad]201 fAreaDefinition = new AreaDefinition(one_ghost_passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
[d7d2da3]202 break;
203 case 3:
[c6321ad]204 fAreaDefinition = new AreaDefinition(passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
[d7d2da3]205 break;
206 case 4:
[c6321ad]207 fAreaDefinition = new AreaDefinition(VoronoiAreaSpec(fEffectiveRfact));
[d7d2da3]208 break;
209 case 5:
[c6321ad]210 fAreaDefinition = new AreaDefinition(active_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
[d7d2da3]211 break;
212 default:
213 case 0:
214 fAreaDefinition = 0;
215 break;
216 }
217
218 switch(fJetAlgorithm)
219 {
[8336b6e]220 case 1:
[c6321ad]221 plugin = new CDFJetCluPlugin(fSeedThreshold, fConeRadius, fAdjacencyCut, fMaxIterations, fIratch, fOverlapThreshold);
222 fDefinition = new JetDefinition(plugin);
[d7d2da3]223 break;
224 case 2:
[c6321ad]225 plugin = new CDFMidPointPlugin(fSeedThreshold, fConeRadius, fConeAreaFraction, fMaxPairSize, fMaxIterations, fOverlapThreshold);
226 fDefinition = new JetDefinition(plugin);
[d7d2da3]227 break;
228 case 3:
[c6321ad]229 plugin = new SISConePlugin(fConeRadius, fOverlapThreshold, fMaxIterations, fJetPTMin);
230 fDefinition = new JetDefinition(plugin);
[d7d2da3]231 break;
232 case 4:
[c6321ad]233 fDefinition = new JetDefinition(kt_algorithm, fParameterR);
[d7d2da3]234 break;
235 case 5:
[c6321ad]236 fDefinition = new JetDefinition(cambridge_algorithm, fParameterR);
[d7d2da3]237 break;
238 default:
239 case 6:
[c6321ad]240 fDefinition = new JetDefinition(antikt_algorithm, fParameterR);
[d7d2da3]241 break;
[9687203]242 case 7:
[c6321ad]243 recomb = new WinnerTakeAllRecombiner();
244 fDefinition = new JetDefinition(antikt_algorithm, fParameterR, recomb, Best);
[9687203]245 break;
246 case 8:
[7278220]247 fNjettinessPlugin = new NjettinessPlugin(fN, Njettiness::wta_kt_axes, Njettiness::unnormalized_cutoff_measure, fBeta, fRcutOff);
248 fDefinition = new JetDefinition(fNjettinessPlugin);
[9687203]249 break;
[f319c1d]250 case 9:
251 fValenciaPlugin = new ValenciaPlugin(fParameterR, fBeta, fGamma);
252 fDefinition = new JetDefinition(fValenciaPlugin);
253 break;
254
[d7d2da3]255 }
[8336b6e]256
[8abab33]257
258
[d7d2da3]259 fPlugin = plugin;
[9687203]260 fRecomb = recomb;
[e4c3fef]261
[d7d2da3]262 ClusterSequence::print_banner();
[8336b6e]263
[7278220]264 if(fComputeRho && fAreaDefinition)
265 {
266 // read eta ranges
267
268 param = GetParam("RhoEtaRange");
269 size = param.GetSize();
270
271 fEstimators.clear();
272 for(i = 0; i < size/2; ++i)
273 {
274 etaMin = param[i*2].GetDouble();
275 etaMax = param[i*2 + 1].GetDouble();
[58dcd2e]276 estimatorStruct.estimator = new JetMedianBackgroundEstimator(SelectorRapRange(etaMin, etaMax), *fDefinition, *fAreaDefinition);
[7278220]277 estimatorStruct.etaMin = etaMin;
278 estimatorStruct.etaMax = etaMax;
279 fEstimators.push_back(estimatorStruct);
280 }
281 }
282
[d7d2da3]283 // import input array
284
285 fInputArray = ImportArray(GetString("InputArray", "Calorimeter/towers"));
286 fItInputArray = fInputArray->MakeIterator();
287
288 // create output arrays
289
290 fOutputArray = ExportArray(GetString("OutputArray", "jets"));
291 fRhoOutputArray = ExportArray(GetString("RhoOutputArray", "rho"));
292}
293
294//------------------------------------------------------------------------------
295
296void FastJetFinder::Finish()
297{
[7278220]298 vector< TEstimatorStruct >::iterator itEstimators;
299
300 for(itEstimators = fEstimators.begin(); itEstimators != fEstimators.end(); ++itEstimators)
301 {
302 if(itEstimators->estimator) delete itEstimators->estimator;
303 }
304
[d7d2da3]305 if(fItInputArray) delete fItInputArray;
306 if(fDefinition) delete fDefinition;
307 if(fAreaDefinition) delete fAreaDefinition;
308 if(fPlugin) delete static_cast<JetDefinition::Plugin*>(fPlugin);
[9687203]309 if(fRecomb) delete static_cast<JetDefinition::Recombiner*>(fRecomb);
310 if(fNjettinessPlugin) delete static_cast<JetDefinition::Plugin*>(fNjettinessPlugin);
[8abab33]311 if(fAxesDef) delete fAxesDef;
312 if(fMeasureDef) delete fMeasureDef;
[f319c1d]313 if(fValenciaPlugin) delete static_cast<JetDefinition::Plugin*>(fValenciaPlugin);
314
[d7d2da3]315}
316
317//------------------------------------------------------------------------------
318
319void FastJetFinder::Process()
320{
321 Candidate *candidate, *constituent;
322 TLorentzVector momentum;
[df04eb1]323
[d7d2da3]324 Double_t deta, dphi, detaMax, dphiMax;
[df04eb1]325 Double_t time, timeWeight;
[ce4feac]326 Int_t number, ncharged, nneutrals;
[5a44a72]327 Int_t charge;
[df04eb1]328 Double_t rho = 0.0;
[d7d2da3]329 PseudoJet jet, area;
330 ClusterSequence *sequence;
[de6d698]331 vector< PseudoJet > inputList, outputList, subjets;
[7278220]332 vector< PseudoJet >::iterator itInputList, itOutputList;
333 vector< TEstimatorStruct >::iterator itEstimators;
[e9c0d73]334 Double_t excl_ymerge23 = 0.0;
335 Double_t excl_ymerge34 = 0.0;
336 Double_t excl_ymerge45 = 0.0;
337 Double_t excl_ymerge56 = 0.0;
338
[d7d2da3]339 DelphesFactory *factory = GetFactory();
340
341 inputList.clear();
342
343 // loop over input objects
344 fItInputArray->Reset();
345 number = 0;
346 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
[8336b6e]347 {
[d7d2da3]348 momentum = candidate->Momentum;
349 jet = PseudoJet(momentum.Px(), momentum.Py(), momentum.Pz(), momentum.E());
350 jet.set_user_index(number);
351 inputList.push_back(jet);
352 ++number;
353 }
[8336b6e]354
355 // construct jets
[d7d2da3]356 if(fAreaDefinition)
357 {
358 sequence = new ClusterSequenceArea(inputList, *fDefinition, *fAreaDefinition);
359 }
360 else
361 {
362 sequence = new ClusterSequence(inputList, *fDefinition);
[8336b6e]363 }
[d7d2da3]364
365 // compute rho and store it
366 if(fComputeRho && fAreaDefinition)
367 {
[7278220]368 for(itEstimators = fEstimators.begin(); itEstimators != fEstimators.end(); ++itEstimators)
[8336b6e]369 {
[7278220]370 itEstimators->estimator->set_particles(inputList);
371 rho = itEstimators->estimator->rho();
[8336b6e]372
373 candidate = factory->NewCandidate();
374 candidate->Momentum.SetPtEtaPhiE(rho, 0.0, 0.0, rho);
[7278220]375 candidate->Edges[0] = itEstimators->etaMin;
376 candidate->Edges[1] = itEstimators->etaMax;
[8336b6e]377 fRhoOutputArray->Add(candidate);
378 }
[d7d2da3]379 }
[8336b6e]380
[d7d2da3]381 outputList.clear();
382
[e9c0d73]383
[f319c1d]384 if(fExclusiveClustering)
385 {
[e9c0d73]386 outputList = sorted_by_pt(sequence->exclusive_jets( fNJets ));
387
388 excl_ymerge23 = sequence->exclusive_ymerge( 2 );
389 excl_ymerge34 = sequence->exclusive_ymerge( 3 );
390 excl_ymerge45 = sequence->exclusive_ymerge( 4 );
391 excl_ymerge56 = sequence->exclusive_ymerge( 5 );
[f319c1d]392 }
393 else
394 {
395 outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin));
396 }
[9687203]397
[d7d2da3]398 // loop over all jets and export them
399 detaMax = 0.0;
400 dphiMax = 0.0;
[ce4feac]401
[d7d2da3]402 for(itOutputList = outputList.begin(); itOutputList != outputList.end(); ++itOutputList)
403 {
[d244bc9]404 jet = *itOutputList;
405 if(fJetAlgorithm == 7) jet = join(jet.constituents());
[df04eb1]406
[d244bc9]407 momentum.SetPxPyPzE(jet.px(), jet.py(), jet.pz(), jet.E());
[df04eb1]408
[d7d2da3]409 area.reset(0.0, 0.0, 0.0, 0.0);
410 if(fAreaDefinition) area = itOutputList->area_4vector();
411
[e9c0d73]412
413
[d7d2da3]414 candidate = factory->NewCandidate();
415
[df04eb1]416 time = 0.0;
417 timeWeight = 0.0;
[22dc7fd]418
[5a44a72]419 charge = 0;
420
[ce4feac]421 ncharged = 0;
422 nneutrals = 0;
423
[d7d2da3]424 inputList.clear();
425 inputList = sequence->constituents(*itOutputList);
[e4c3fef]426
[d7d2da3]427 for(itInputList = inputList.begin(); itInputList != inputList.end(); ++itInputList)
428 {
[29f8a06]429 if(itInputList->user_index() < 0) continue;
[d7d2da3]430 constituent = static_cast<Candidate*>(fInputArray->At(itInputList->user_index()));
431
432 deta = TMath::Abs(momentum.Eta() - constituent->Momentum.Eta());
433 dphi = TMath::Abs(momentum.DeltaPhi(constituent->Momentum));
434 if(deta > detaMax) detaMax = deta;
435 if(dphi > dphiMax) dphiMax = dphi;
[e4c3fef]436
[ce4feac]437 if(constituent->Charge == 0) nneutrals++;
438 else ncharged++;
439
[22dc7fd]440 time += TMath::Sqrt(constituent->Momentum.E())*(constituent->Position.T());
[df04eb1]441 timeWeight += TMath::Sqrt(constituent->Momentum.E());
[e4c3fef]442
[5a44a72]443 charge += constituent->Charge;
444
[d7d2da3]445 candidate->AddCandidate(constituent);
446 }
[e4c3fef]447
[d7d2da3]448 candidate->Momentum = momentum;
[df04eb1]449 candidate->Position.SetT(time/timeWeight);
[d7d2da3]450 candidate->Area.SetPxPyPzE(area.px(), area.py(), area.pz(), area.E());
451
452 candidate->DeltaEta = detaMax;
453 candidate->DeltaPhi = dphiMax;
[5a44a72]454 candidate->Charge = charge;
[ce4feac]455 candidate->NNeutrals = nneutrals;
456 candidate->NCharged = ncharged;
[e9c0d73]457
458
459 //for exclusive clustering, access y_n,n+1 as exclusive_ymerge (fNJets);
460 candidate->ExclYmerge23 = excl_ymerge23;
461 candidate->ExclYmerge34 = excl_ymerge34;
462 candidate->ExclYmerge45 = excl_ymerge45;
463 candidate->ExclYmerge56 = excl_ymerge56;
[ce4feac]464
[de6d698]465 //------------------------------------
466 // Trimming
467 //------------------------------------
[d7d2da3]468
[8abab33]469 if(fComputeTrimming)
[de6d698]470 {
471
472 fastjet::Filter trimmer(fastjet::JetDefinition(fastjet::kt_algorithm,fRTrim),fastjet::SelectorPtFractionMin(fPtFracTrim));
473 fastjet::PseudoJet trimmed_jet = trimmer(*itOutputList);
474
475 trimmed_jet = join(trimmed_jet.constituents());
476
477 candidate->TrimmedP4[0].SetPtEtaPhiM(trimmed_jet.pt(), trimmed_jet.eta(), trimmed_jet.phi(), trimmed_jet.m());
478
479 // four hardest subjets
480 subjets.clear();
481 subjets = trimmed_jet.pieces();
482 subjets = sorted_by_pt(subjets);
483
484 candidate->NSubJetsTrimmed = subjets.size();
485
[8abab33]486 for (size_t i = 0; i < subjets.size() and i < 4; i++)
487 {
488 if(subjets.at(i).pt() < 0) continue ;
489 candidate->TrimmedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m());
[de6d698]490 }
491 }
492
493
494 //------------------------------------
495 // Pruning
496 //------------------------------------
497
498
499 if(fComputePruning)
500 {
501
502 fastjet::Pruner pruner(fastjet::JetDefinition(fastjet::cambridge_algorithm,fRPrun),fZcutPrun,fRcutPrun);
503 fastjet::PseudoJet pruned_jet = pruner(*itOutputList);
504
505 candidate->PrunedP4[0].SetPtEtaPhiM(pruned_jet.pt(), pruned_jet.eta(), pruned_jet.phi(), pruned_jet.m());
506
507 // four hardest subjet
508 subjets.clear();
509 subjets = pruned_jet.pieces();
510 subjets = sorted_by_pt(subjets);
511
512 candidate->NSubJetsPruned = subjets.size();
513
[8abab33]514 for (size_t i = 0; i < subjets.size() and i < 4; i++)
515 {
516 if(subjets.at(i).pt() < 0) continue ;
517 candidate->PrunedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m());
[de6d698]518 }
519
520 }
521
522 //------------------------------------
523 // SoftDrop
524 //------------------------------------
525
526 if(fComputeSoftDrop)
527 {
528
[8abab33]529 contrib::SoftDrop softDrop(fBetaSoftDrop,fSymmetryCutSoftDrop,fR0SoftDrop);
[de6d698]530 fastjet::PseudoJet softdrop_jet = softDrop(*itOutputList);
531
532 candidate->SoftDroppedP4[0].SetPtEtaPhiM(softdrop_jet.pt(), softdrop_jet.eta(), softdrop_jet.phi(), softdrop_jet.m());
533
534 // four hardest subjet
535
536 subjets.clear();
537 subjets = softdrop_jet.pieces();
538 subjets = sorted_by_pt(subjets);
539 candidate->NSubJetsSoftDropped = softdrop_jet.pieces().size();
540
[ba75867]541 candidate->SoftDroppedJet = candidate->SoftDroppedP4[0];
542
[8abab33]543 for (size_t i = 0; i < subjets.size() and i < 4; i++)
544 {
545 if(subjets.at(i).pt() < 0) continue ;
546 candidate->SoftDroppedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m());
[ba75867]547 if(i==0) candidate->SoftDroppedSubJet1 = candidate->SoftDroppedP4[i+1];
548 if(i==1) candidate->SoftDroppedSubJet2 = candidate->SoftDroppedP4[i+1];
[de6d698]549 }
550 }
551
[9687203]552 // --- compute N-subjettiness with N = 1,2,3,4,5 ----
[e4c3fef]553
[9687203]554 if(fComputeNsubjettiness)
555 {
[8abab33]556
557 Nsubjettiness nSub1(1, *fAxesDef, *fMeasureDef);
558 Nsubjettiness nSub2(2, *fAxesDef, *fMeasureDef);
559 Nsubjettiness nSub3(3, *fAxesDef, *fMeasureDef);
560 Nsubjettiness nSub4(4, *fAxesDef, *fMeasureDef);
561 Nsubjettiness nSub5(5, *fAxesDef, *fMeasureDef);
562
[e4c3fef]563 candidate->Tau[0] = nSub1(*itOutputList);
564 candidate->Tau[1] = nSub2(*itOutputList);
565 candidate->Tau[2] = nSub3(*itOutputList);
566 candidate->Tau[3] = nSub4(*itOutputList);
567 candidate->Tau[4] = nSub5(*itOutputList);
[8abab33]568
[9687203]569 }
570
[d7d2da3]571 fOutputArray->Add(candidate);
572 }
573 delete sequence;
574}
Note: See TracBrowser for help on using the repository browser.