Fork me on GitHub

source: git/external/fastjet/contribs/Nsubjettiness/NjettinessPlugin.cc@ 51cabe8

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

upgrade FastJet to version 3.1.0-beta.1, upgrade Nsubjettiness to version 2.1.0, add SoftKiller version 1.0.0

  • Property mode set to 100644
File size: 3.3 KB
Line 
1// Nsubjettiness Package
2// Questions/Comments? jthaler@jthaler.net
3//
4// Copyright (c) 2011-14
5// Jesse Thaler, Ken Van Tilburg, Christopher K. Vermilion, and TJ Wilkason
6//
7// $Id: NjettinessPlugin.cc 663 2014-06-03 21:26:41Z jthaler $
8//----------------------------------------------------------------------
9// This file is part of FastJet contrib.
10//
11// It is free software; you can redistribute it and/or modify it under
12// the terms of the GNU General Public License as published by the
13// Free Software Foundation; either version 2 of the License, or (at
14// your option) any later version.
15//
16// It is distributed in the hope that it will be useful, but WITHOUT
17// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
18// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
19// License for more details.
20//
21// You should have received a copy of the GNU General Public License
22// along with this code. If not, see <http://www.gnu.org/licenses/>.
23//----------------------------------------------------------------------
24
25#include "NjettinessPlugin.hh"
26
27FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
28
29namespace contrib{
30
31
32
33std::string NjettinessPlugin::description() const {return "N-jettiness jet finder";}
34
35
36// Clusters the particles according to the Njettiness jet algorithm
37// Apologies for the complication with this code, but we need to make
38// a fake jet clustering tree. The partitioning is done by getPartitionList
39void NjettinessPlugin::run_clustering(ClusterSequence& cs) const
40{
41 std::vector<fastjet::PseudoJet> particles = cs.jets();
42
43 // HACK: remove area information from particles (in case this is called by
44 // a ClusterSequenceArea. Will be fixed in a future FastJet release)
45 for (unsigned i = 0; i < particles.size(); i++) {
46 particles[i].set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>());
47 }
48
49
50 _njettinessFinder.getTau(_N, particles);
51
52 std::vector<std::list<int> > partition = _njettinessFinder.getPartitionList(particles);
53
54 std::vector<fastjet::PseudoJet> jet_indices_for_extras;
55
56 // output clusterings for each jet
57 for (size_t i0 = 0; i0 < partition.size(); ++i0) {
58 size_t i = partition.size() - 1 - i0; // reversed order of reading to match axes order
59 std::list<int>& indices = partition[i];
60 if (indices.size() == 0) continue;
61 while (indices.size() > 1) {
62 int merge_i = indices.back(); indices.pop_back();
63 int merge_j = indices.back(); indices.pop_back();
64 int newIndex;
65 double fakeDij = -1.0;
66
67 cs.plugin_record_ij_recombination(merge_i, merge_j, fakeDij, newIndex);
68
69 indices.push_back(newIndex);
70 }
71 double fakeDib = -1.0;
72
73 int finalJet = indices.back();
74 cs.plugin_record_iB_recombination(finalJet, fakeDib);
75 jet_indices_for_extras.push_back(cs.jets()[finalJet]); // Get the four vector for the final jets to compare later.
76 }
77
78 //HACK: Re-reverse order of reading to match CS order
79 reverse(jet_indices_for_extras.begin(),jet_indices_for_extras.end());
80
81 NjettinessExtras * extras = new NjettinessExtras(_njettinessFinder.currentTauComponents(),jet_indices_for_extras,_njettinessFinder.currentAxes());
82 cs.plugin_associate_extras(std::auto_ptr<ClusterSequence::Extras>(extras));
83
84}
85
86
87} // namespace contrib
88
89FASTJET_END_NAMESPACE
Note: See TracBrowser for help on using the repository browser.