FastJet  3.1.2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
ATLASConePlugin.hh
1 //FJSTARTHEADER
2 // $Id: ATLASConePlugin.hh 3433 2014-07-23 08:17:03Z salam $
3 //
4 // Copyright (c) 2007-2014, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5 //
6 //----------------------------------------------------------------------
7 // This file is part of FastJet.
8 //
9 // FastJet is free software; you can redistribute it and/or modify
10 // it under the terms of the GNU General Public License as published by
11 // the Free Software Foundation; either version 2 of the License, or
12 // (at your option) any later version.
13 //
14 // The algorithms that underlie FastJet have required considerable
15 // development. They are described in the original FastJet paper,
16 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17 // FastJet as part of work towards a scientific publication, please
18 // quote the version you use and include a citation to the manual and
19 // optionally also to hep-ph/0512210.
20 //
21 // FastJet is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28 //----------------------------------------------------------------------
29 //FJENDHEADER
30 
31 // Note on the implementation:
32 // this implementation of the ATLAS Cone is based on the SpartyJet
33 // v2.20.0 implementation. See README for details
34 
35 #ifndef __ATLASCONEPLUGIN_HH__
36 #define __ATLASCONEPLUGIN_HH__
37 
38 #include "fastjet/JetDefinition.hh"
39 
40 // questionable whether this should be in fastjet namespace or not...
41 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
42 
43 // forward declaration to reduce includes
44 class PseudoJet;
45 
46 //----------------------------------------------------------------------
47 //
48 /// @ingroup plugins
49 /// \class ATLASConePlugin
50 /// Implementation of the ATLAS Cone (plugin for fastjet v2.4 upwards)
52 public:
53  /// Main constructor for the ATLASCone Plugin class.
54  ///
55  /// Apparently the default parameters in the ATLAS software are the
56  /// ones used here. SpartyJet uses a radius of 0.7, a seed threshold
57  /// of 1 GeV and an overlap threshold of 0.75
58  /// For the ATLAS SW defaults, see
59  /// http://atlas-sw.cern.ch/cgi-bin/viewcvs-atlas.cgi/groups/JetRoutines/SpartyJet/atlas/
60  /// in the JetdoneFinderTools.cxx (rev1.1) and JetSplitMergeTool.cxx (rev1.1)
61  /// For SpartyJet, see atlas/ConeFinderTool.h
62  ///
63  /// Finally, to agree with FastJet standards, we do not specify a default R,
64  /// that in the ATLAS code is 0.7
65  ATLASConePlugin (double radius, double seedPt_in=2.0, double f_in=0.5)
66  : _radius(radius), _seedPt(seedPt_in), _f(f_in){}
67 
68  /// copy constructor
69  ATLASConePlugin (const ATLASConePlugin & plugin) {
70  *this = plugin;
71  }
72 
73  // the things that are required by base class
74  virtual std::string description () const;
75  virtual void run_clustering(ClusterSequence &) const;
76 
77  /// the plugin mechanism's standard way of accessing the jet radius
78  /// here we return the R of the last alg in the list
79  virtual double R() const {return _radius;}
80 
81  // access to the other parameters
82  /// seed threshold
83  double seedPt() const {return _seedPt;}
84 
85  /// split-merge overlap threshold
86  double f() const {return _f;}
87 
88 private:
89 
90  double _radius; ///< the cone radius
91  double _seedPt; ///< the pt seed threshold used in stable-cone search
92  double _f; ///< the overlap thresholod used in the split-merge
93 
94  static bool _first_time;
95 
96  /// print a banner for reference to the 3rd-party code
97  void _print_banner(std::ostream *ostr) const;
98 };
99 
100 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
101 
102 #endif // __ATLASCONEPLUGIN_HH__
103 
deals with clustering
virtual double R() const
the plugin mechanism's standard way of accessing the jet radius here we return the R of the last alg ...
double seedPt() const
seed threshold
ATLASConePlugin(double radius, double seedPt_in=2.0, double f_in=0.5)
Main constructor for the ATLASCone Plugin class.
a class that allows a user to introduce their own "plugin" jet finder
double f() const
split-merge overlap threshold
Implementation of the ATLAS Cone (plugin for fastjet v2.4 upwards)
ATLASConePlugin(const ATLASConePlugin &plugin)
copy constructor