1 | //STARTHEADER |
---|
2 | // $Id: ATLASConePlugin.hh 2758 2011-11-24 08:31:58Z soyez $ |
---|
3 | // |
---|
4 | // Copyright (c) 2007-2011, 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 and are described in hep-ph/0512210. If you use |
---|
16 | // FastJet as part of work towards a scientific publication, please |
---|
17 | // include a citation to the FastJet paper. |
---|
18 | // |
---|
19 | // FastJet is distributed in the hope that it will be useful, |
---|
20 | // but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
21 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
22 | // GNU General Public License for more details. |
---|
23 | // |
---|
24 | // You should have received a copy of the GNU General Public License |
---|
25 | // along with FastJet. If not, see <http://www.gnu.org/licenses/>. |
---|
26 | //---------------------------------------------------------------------- |
---|
27 | //ENDHEADER |
---|
28 | |
---|
29 | // Note on the implementation: |
---|
30 | // this implementation of the ATLAS Cone is based on the SpartyJet |
---|
31 | // v2.20.0 implementation. See README for details |
---|
32 | |
---|
33 | #ifndef __ATLASCONEPLUGIN_HH__ |
---|
34 | #define __ATLASCONEPLUGIN_HH__ |
---|
35 | |
---|
36 | #include "fastjet/JetDefinition.hh" |
---|
37 | |
---|
38 | // questionable whether this should be in fastjet namespace or not... |
---|
39 | FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh |
---|
40 | |
---|
41 | // forward declaration to reduce includes |
---|
42 | class PseudoJet; |
---|
43 | |
---|
44 | //---------------------------------------------------------------------- |
---|
45 | // |
---|
46 | /// @ingroup plugins |
---|
47 | /// \class ATLASConePlugin |
---|
48 | /// Implementation of the ATLAS Cone (plugin for fastjet v2.4 upwards) |
---|
49 | class ATLASConePlugin : public JetDefinition::Plugin { |
---|
50 | public: |
---|
51 | /// Main constructor for the ATLASCone Plugin class. |
---|
52 | /// |
---|
53 | /// Apparently the default parameters in the ATLAS software are the |
---|
54 | /// ones used here. SpartyJet uses a radius of 0.7, a seed threshold |
---|
55 | /// of 1 GeV and an overlap threshold of 0.75 |
---|
56 | /// For the ATLAS SW defaults, see |
---|
57 | /// http://atlas-sw.cern.ch/cgi-bin/viewcvs-atlas.cgi/groups/JetRoutines/SpartyJet/atlas/ |
---|
58 | /// in the JetdoneFinderTools.cxx (rev1.1) and JetSplitMergeTool.cxx (rev1.1) |
---|
59 | /// For SpartyJet, see atlas/ConeFinderTool.h |
---|
60 | /// |
---|
61 | /// Finally, to agree with FastJet standards, we do not specify a default R, |
---|
62 | /// that in the ATLAS code is 0.7 |
---|
63 | ATLASConePlugin (double radius, double seedPt_in=2.0, double f_in=0.5) |
---|
64 | : _radius(radius), _seedPt(seedPt_in), _f(f_in){} |
---|
65 | |
---|
66 | /// copy constructor |
---|
67 | ATLASConePlugin (const ATLASConePlugin & plugin) { |
---|
68 | *this = plugin; |
---|
69 | } |
---|
70 | |
---|
71 | // the things that are required by base class |
---|
72 | virtual std::string description () const; |
---|
73 | virtual void run_clustering(ClusterSequence &) const; |
---|
74 | |
---|
75 | /// the plugin mechanism's standard way of accessing the jet radius |
---|
76 | /// here we return the R of the last alg in the list |
---|
77 | virtual double R() const {return _radius;} |
---|
78 | |
---|
79 | // access to the other parameters |
---|
80 | /// seed threshold |
---|
81 | double seedPt() const {return _seedPt;} |
---|
82 | |
---|
83 | /// split-merge overlap threshold |
---|
84 | double f() const {return _f;} |
---|
85 | |
---|
86 | private: |
---|
87 | |
---|
88 | double _radius; ///< the cone radius |
---|
89 | double _seedPt; ///< the pt seed threshold used in stable-cone search |
---|
90 | double _f; ///< the overlap thresholod used in the split-merge |
---|
91 | |
---|
92 | static bool _first_time; |
---|
93 | |
---|
94 | /// print a banner for reference to the 3rd-party code |
---|
95 | void _print_banner(std::ostream *ostr) const; |
---|
96 | }; |
---|
97 | |
---|
98 | FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh |
---|
99 | |
---|
100 | #endif // __ATLASCONEPLUGIN_HH__ |
---|
101 | |
---|