Fork me on GitHub

source: git/external/fastjet/contribs/Nsubjettiness/AxesFinder.hh@ 50d7259

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 50d7259 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: 13.6 KB
RevLine 
[9687203]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//
[35cdc46]7// $Id: AxesFinder.hh 678 2014-06-12 20:43:03Z jthaler $
[9687203]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
26#ifndef __FASTJET_CONTRIB_AXESFINDER_HH__
27#define __FASTJET_CONTRIB_AXESFINDER_HH__
28
29#include "WinnerTakeAllRecombiner.hh"
30#include "MeasureFunction.hh"
31
32#include "fastjet/PseudoJet.hh"
33#include "fastjet/ClusterSequence.hh"
34#include "fastjet/JetDefinition.hh"
35
36#include <cmath>
37#include <vector>
38#include <list>
39
40FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41
42namespace contrib{
43
44///////
45//
46// Axes Finder Options
47//
48///////
49
50//------------------------------------------------------------------------
51/// \class AxesFinder
52// This is the base class for all axes finders. These axes are used along with the MeasureFunctions to calculate
53// tau_N. There are different implementations of axes finding that are defined in derived classes below.
54class AxesFinder {
55
56public:
57
[35cdc46]58 // This function should be overloaded, and updates the seedAxes to return new axes
59 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets,
60 const std::vector<fastjet::PseudoJet>& inputs,
61 const std::vector<fastjet::PseudoJet>& seedAxes) const = 0;
62 // convenient shorthand for squaring
63 static inline double sq(double x) {return x*x;}
64
65 //virtual destructor
66 virtual ~AxesFinder(){}
[9687203]67
68};
69
[35cdc46]70
[9687203]71//------------------------------------------------------------------------
72/// \class AxesFinderFromExclusiveJetDefinition
73// This class finds axes by clustering the particles and then finding the exclusive jets. This can be implemented
74// with different jet algorithms.
75class AxesFinderFromExclusiveJetDefinition : public AxesFinder {
76
[35cdc46]77public:
78 AxesFinderFromExclusiveJetDefinition(fastjet::JetDefinition def)
79 : _def(def) {}
80
81 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets,
82 const std::vector <fastjet::PseudoJet> & inputs,
83 const std::vector<fastjet::PseudoJet>& /*seedAxes*/) const {
84 fastjet::ClusterSequence jet_clust_seq(inputs, _def);
85 return jet_clust_seq.exclusive_jets(n_jets);
86 }
87
88private:
89 fastjet::JetDefinition _def;
90
[9687203]91};
92
93//------------------------------------------------------------------------
94/// \class AxesFinderFromWTA_KT
95// This class finds axes by finding the exlusive jets after clustering according to a kT algorithm and a
96// winner take all recombination scheme.
97class AxesFinderFromWTA_KT : public AxesFinderFromExclusiveJetDefinition {
[35cdc46]98
99public:
100 AxesFinderFromWTA_KT()
101 : AxesFinderFromExclusiveJetDefinition(
102 fastjet::JetDefinition(fastjet::kt_algorithm,
103 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
104 &_recomb,
105 fastjet::Best)) {}
106
107private:
108 const WinnerTakeAllRecombiner _recomb;
109
110};
[9687203]111
112//------------------------------------------------------------------------
113/// \class AxesFinderFromWTA_CA
114// This class finds axes by finding the exlusive jets after clustering according to a CA algorithm and a
115// winner take all recombination scheme.
116class AxesFinderFromWTA_CA : public AxesFinderFromExclusiveJetDefinition {
[35cdc46]117public:
118 AxesFinderFromWTA_CA()
119 : AxesFinderFromExclusiveJetDefinition(
120 fastjet::JetDefinition(fastjet::cambridge_algorithm,
121 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
122 &_recomb,
123 fastjet::Best)) {}
124
125private:
126 const WinnerTakeAllRecombiner _recomb;
[9687203]127};
128
129
130//------------------------------------------------------------------------
131/// \class AxesFinderFromKT
132// This class finds axes by finding the exlusive jets after clustering according to a kT algorithm and a
133// E_scheme recombination.
134class AxesFinderFromKT : public AxesFinderFromExclusiveJetDefinition {
[35cdc46]135public:
136 AxesFinderFromKT()
137 : AxesFinderFromExclusiveJetDefinition(
138 fastjet::JetDefinition(fastjet::kt_algorithm,
139 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
140 fastjet::E_scheme,
141 fastjet::Best)) {}
[9687203]142};
143
144//------------------------------------------------------------------------
145/// \class AxesFinderFromCA
146// This class finds axes by finding the exlusive jets after clustering according to a CA algorithm and a
147// E_scheme recombination.
148class AxesFinderFromCA : public AxesFinderFromExclusiveJetDefinition {
[35cdc46]149public:
150 AxesFinderFromCA()
151 : AxesFinderFromExclusiveJetDefinition(
152 fastjet::JetDefinition(fastjet::cambridge_algorithm,
153 fastjet::JetDefinition::max_allowable_R, //maximum jet radius constant
154 fastjet::E_scheme,
155 fastjet::Best)) {}
[9687203]156};
157
158
159//------------------------------------------------------------------------
160/// \class AxesFinderFromHardestJetDefinition
161// This class finds axes by clustering the particles and then finding the n hardest inclusive jets.
162// This can be implemented with different jet algorithms.
163class AxesFinderFromHardestJetDefinition : public AxesFinder {
[35cdc46]164public:
165 AxesFinderFromHardestJetDefinition(fastjet::JetDefinition def)
166 : _def(def) {}
[9687203]167
[35cdc46]168 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets,
169 const std::vector <fastjet::PseudoJet> & inputs,
170 const std::vector<fastjet::PseudoJet>& /*seedAxes*/) const {
171 fastjet::ClusterSequence jet_clust_seq(inputs, _def);
172 std::vector<fastjet::PseudoJet> myJets = sorted_by_pt(jet_clust_seq.inclusive_jets());
173 myJets.resize(n_jets); // only keep n hardest
174 return myJets;
175 }
176
177private:
178 fastjet::JetDefinition _def;
[9687203]179};
180
181//------------------------------------------------------------------------
182/// \class AxesFinderFromAntiKT
183// This class finds axes by finding the n hardest jets after clustering the particles according
184// to an anti kT algorithm and E_scheme.
185class AxesFinderFromAntiKT : public AxesFinderFromHardestJetDefinition {
[35cdc46]186public:
187 AxesFinderFromAntiKT(double R0)
188 : AxesFinderFromHardestJetDefinition(
189 fastjet::JetDefinition(fastjet::antikt_algorithm,
190 R0,fastjet::E_scheme,fastjet::Best)) {}
[9687203]191};
192
193
194//------------------------------------------------------------------------
195/// \class AxesFinderFromUserInput
196// This class allows the user to manually define the axes.
197class AxesFinderFromUserInput : public AxesFinder {
198
[35cdc46]199public:
200 AxesFinderFromUserInput() {}
201
202 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets, const std::vector <fastjet::PseudoJet> & /*inputs*/, const std::vector<fastjet::PseudoJet>& currentAxes) const {
203 assert(currentAxes.size() == (unsigned int) n_jets);
204 (void)(n_jets); // adding this line to fix unused-parameter warning
205 return currentAxes;
206 }
[9687203]207};
208
209//This is a helper class for the Minimum Axes Finders. It is defined later.
210class LightLikeAxis;
211
212
213//------------------------------------------------------------------------
214/// \class AxesFinderFromOnePassMinimization
215// This class defines an AxesFinder that uses Kmeans minimization, but only on a single pass.
216class AxesFinderFromOnePassMinimization : public AxesFinder {
217
[35cdc46]218public:
219
220 // From a startingFinder, try to minimize the unnormalized_measure
221 AxesFinderFromOnePassMinimization(double beta, double Rcutoff)
222 : _precision(0.0001), //hard coded for now
223 _halt(1000), //hard coded for now
224 _beta(beta),
225 _Rcutoff(Rcutoff),
226 _measureFunction(beta, Rcutoff)
227 {}
228
229 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets,
230 const std::vector <fastjet::PseudoJet> & inputJets,
231 const std::vector<fastjet::PseudoJet>& currentAxes) const;
[9687203]232
[35cdc46]233private:
234 double _precision; // Desired precision in axes alignment
235 int _halt; // maximum number of steps per iteration
[9687203]236
[35cdc46]237 double _beta;
238 double _Rcutoff;
[9687203]239
[35cdc46]240 DefaultUnnormalizedMeasureFunction _measureFunction;
241
242 template <int N> std::vector<LightLikeAxis> UpdateAxesFast(const std::vector <LightLikeAxis> & old_axes,
243 const std::vector <fastjet::PseudoJet> & inputJets) const;
244
245 std::vector<LightLikeAxis> UpdateAxes(const std::vector <LightLikeAxis> & old_axes,
246 const std::vector <fastjet::PseudoJet> & inputJets) const;
[9687203]247
248};
249
250
251//------------------------------------------------------------------------
252/// \class AxesFinderFromKmeansMinimization
253// This class finds finds axes by using Kmeans clustering to minimizaiton N-jettiness. Given a first set of
254// starting axes, it updates n times to get as close to the global minimum as possible. This class calls OnePass many times,
255// added noise to the axes.
256class AxesFinderFromKmeansMinimization : public AxesFinder{
257
[35cdc46]258public:
259 AxesFinderFromKmeansMinimization(double beta, double Rcutoff, int n_iterations)
260 : _n_iterations(n_iterations),
261 _noise_range(1.0), // hard coded for the time being
262 _measureFunction(beta, Rcutoff),
263 _onePassFinder(beta, Rcutoff)
264 {}
265
266 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets, const std::vector <fastjet::PseudoJet> & inputJets, const std::vector<fastjet::PseudoJet>& currentAxes) const;
[9687203]267
[35cdc46]268private:
269 int _n_iterations; // Number of iterations to run (0 for no minimization, 1 for one-pass, >>1 for global minimum)
270 double _noise_range; // noise range for random initialization
[9687203]271
[35cdc46]272 DefaultUnnormalizedMeasureFunction _measureFunction; //function to test whether minimum is reached
[9687203]273
[35cdc46]274 AxesFinderFromOnePassMinimization _onePassFinder; //one pass finder that is repeatedly called
275
276 PseudoJet jiggle(const PseudoJet& axis) const;
[9687203]277};
278
279//------------------------------------------------------------------------
280/// \class AxesFinderFromGeometricMinimization
281// This class finds axes by minimizing the Lorentz dot product distance between axes and particles. Given a first set of starting axes,
282// it essentially does stable cone finxing.
283class AxesFinderFromGeometricMinimization : public AxesFinder {
284
[35cdc46]285public:
286 AxesFinderFromGeometricMinimization(double beta, double Rcutoff)
287 : _nAttempts(100),
288 _accuracy(0.000000001),
289 _function(beta,Rcutoff)
290 {
291 if (beta != 2.0) {
292 throw Error("Geometric minimization is currently only defined for beta = 2.0.");
[9687203]293 }
[35cdc46]294 }
295
296 virtual std::vector<fastjet::PseudoJet> getAxes(int n_jets, const std::vector <fastjet::PseudoJet> & particles, const std::vector<fastjet::PseudoJet>& currentAxes) const;
297
298private:
299 double _nAttempts;
300 double _accuracy;
301 GeometricMeasureFunction _function;
302
[9687203]303
304};
305
306//------------------------------------------------------------------------
307/// \class LightLikeAxis
308// This is a helper class for the minimum Axes Finders classes above. It creates a convenient way of defining axes
309// in order to better facilitate calculations.
310class LightLikeAxis {
[35cdc46]311
[9687203]312public:
313 LightLikeAxis() : _rap(0.0), _phi(0.0), _weight(0.0), _mom(0.0) {}
314 LightLikeAxis(double my_rap, double my_phi, double my_weight, double my_mom) :
315 _rap(my_rap), _phi(my_phi), _weight(my_weight), _mom(my_mom) {}
316 double rap() const {return _rap;}
317 double phi() const {return _phi;}
318 double weight() const {return _weight;}
319 double mom() const {return _mom;}
320 void set_rap(double my_set_rap) {_rap = my_set_rap;}
321 void set_phi(double my_set_phi) {_phi = my_set_phi;}
322 void set_weight(double my_set_weight) {_weight = my_set_weight;}
323 void set_mom(double my_set_mom) {_mom = my_set_mom;}
324 void reset(double my_rap, double my_phi, double my_weight, double my_mom) {_rap=my_rap; _phi=my_phi; _weight=my_weight; _mom=my_mom;}
325
326 // return PseudoJet with information
327 fastjet::PseudoJet ConvertToPseudoJet();
328
329 double DistanceSq(const fastjet::PseudoJet& input) const {
330 return DistanceSq(input.rap(),input.phi());
331 }
332
333 double Distance(const fastjet::PseudoJet& input) const {
334 return std::sqrt(DistanceSq(input));
335 }
336
337 double DistanceSq(const LightLikeAxis& input) const {
338 return DistanceSq(input.rap(),input.phi());
339 }
340
341 double Distance(const LightLikeAxis& input) const {
342 return std::sqrt(DistanceSq(input));
343 }
344
[35cdc46]345private:
346 double _rap, _phi, _weight, _mom;
347
348 double DistanceSq(double rap2, double phi2) const {
349 double rap1 = _rap;
350 double phi1 = _phi;
351
352 double distRap = rap1-rap2;
353 double distPhi = std::fabs(phi1-phi2);
354 if (distPhi > M_PI) {distPhi = 2.0*M_PI - distPhi;}
355 return distRap*distRap + distPhi*distPhi;
356 }
357
358 double Distance(double rap2, double phi2) const {
359 return std::sqrt(DistanceSq(rap2,phi2));
360 }
361
362
[9687203]363};
364
365} //namespace contrib
366
367FASTJET_END_NAMESPACE
368
369#endif // __FASTJET_CONTRIB_AXESFINDER_HH__
Note: See TracBrowser for help on using the repository browser.