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