Fork me on GitHub

source: git/external/fastjet/contribs/Nsubjettiness/NjettinessPlugin.cc@ 5babc8b

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5babc8b was 9687203, checked in by mselvaggi <mselvaggi@…>, 11 years ago

added nsubjettiness

  • Property mode set to 100644
File size: 3.2 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//----------------------------------------------------------------------
8// This file is part of FastJet contrib.
9//
10// It is free software; you can redistribute it and/or modify it under
11// the terms of the GNU General Public License as published by the
12// Free Software Foundation; either version 2 of the License, or (at
13// your option) any later version.
14//
15// It is distributed in the hope that it will be useful, but WITHOUT
16// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
17// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
18// License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with this code. If not, see <http://www.gnu.org/licenses/>.
22//----------------------------------------------------------------------
23
24#include "NjettinessPlugin.hh"
25
26FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
27
28namespace contrib{
29
30
31// Constructor with same arguments as Nsubjettiness.
32NjettinessPlugin::NjettinessPlugin(int N, Njettiness::AxesMode axes_mode, Njettiness::MeasureMode measure_mode, double para1, double para2, double para3, double para4)
33 : _N(N), _njettinessFinder(axes_mode, measure_mode, para1, para2, para3, para4) {}
34
35// Old constructor for compatibility
36NjettinessPlugin::NjettinessPlugin(int N, Njettiness::AxesMode mode, double beta, double R0, double Rcutoff)
37 : _N(N), _njettinessFinder(mode, Njettiness::normalized_cutoff_measure, beta, R0, Rcutoff) {}
38
39std::string NjettinessPlugin::description() const {return "N-jettiness jet finder";}
40
41// Clusters the particles according to the Njettiness jet algorithm
42// TODO: this code should be revisited to see if if can be made more clear.
43void NjettinessPlugin::run_clustering(ClusterSequence& cs) const
44{
45 std::vector<fastjet::PseudoJet> particles = cs.jets();
46 _njettinessFinder.getTau(_N, particles);
47 std::vector<std::list<int> > partition = _njettinessFinder.getPartition(particles);
48
49 std::vector<fastjet::PseudoJet> jet_indices_for_extras;
50
51 // output clusterings for each jet
52 for (size_t i = 0; i < partition.size(); ++i) {
53 std::list<int>& indices = partition[i];
54 if (indices.size() == 0) continue;
55 while (indices.size() > 1) {
56 int merge_i = indices.back(); indices.pop_back();
57 int merge_j = indices.back(); indices.pop_back();
58 int newIndex;
59 double fakeDij = -1.0;
60
61 cs.plugin_record_ij_recombination(merge_i, merge_j, fakeDij, newIndex);
62
63 indices.push_back(newIndex);
64 }
65 double fakeDib = -1.0;
66
67 int finalJet = indices.back();
68 cs.plugin_record_iB_recombination(finalJet, fakeDib);
69 jet_indices_for_extras.push_back(cs.jets()[finalJet]); // Get the four vector for the final jets to compare later.
70 }
71
72 NjettinessExtras * extras = new NjettinessExtras(_njettinessFinder.currentTauComponents(),jet_indices_for_extras,_njettinessFinder.currentAxes());
73 cs.plugin_associate_extras(std::auto_ptr<ClusterSequence::Extras>(extras));
74
75}
76
77
78} // namespace contrib
79
80FASTJET_END_NAMESPACE
Note: See TracBrowser for help on using the repository browser.