FastJet 3.4.1
PxConePlugin.hh
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2023, 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#ifndef __PXCONEPLUGIN_HH__
32#define __PXCONEPLUGIN_HH__
33
34#include "fastjet/JetDefinition.hh"
35#include "fastjet/internal/thread_safety_helpers.hh" // helpers to write transparent code w&wo C++11 features
36
37// questionable whether this should be in fastjet namespace or not...
38
39FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
40
41//----------------------------------------------------------------------
42//
43/// @ingroup plugins
44/// \class PxConePlugin
45/// Implementation of the PxCone algorithm (plugin for fastjet v2.1 upwards)
46///
47/// PxConePlugin is a plugin for fastjet (v2.1 upwards) that provides
48/// an interface to the fortran pxcone iterative cone algorithm with
49/// midpoint seeds.
50///
51/// Pxcone was written by Luis del Pozo and Michael H. Seymour. It is
52/// not a "supported" program, so if you encounter problems, you are
53/// on your own...
54///
55/// Note that pxcone sometimes encounters non-stable iterations; in
56/// such cases it returns an error -- the plugin propagates this by
57/// throwing a fastjet::Error exception; if the user wishes to have
58/// robust code, they should catch this exception.
59///
60/// Pxcone has a hard-coded limit (by default 4000) on the maximum
61/// number of particles and protojets; if the number of particles or
62/// protojets exceeds this, again a fastjet::Error exception will be
63/// thrown.
64///
65/// The functionality of pxcone is described at
66/// http://www.hep.man.ac.uk/u/wplano/ConeJet.ps
67//
68//----------------------------------------------------------------------
70public:
71
72 /// constructor for the PxConePlugin, whose arguments have the
73 /// following meaning:
74 ///
75 /// - the cone_radius is as usual in cone algorithms
76 ///
77 /// - stables cones (protojets) below min_jet_energy are discarded
78 /// before calling the splitting procedure to resolve overlaps
79 /// (called epslon in pxcone).
80 ///
81 /// - when two protojets overlap, if
82 /// (overlapping_Et)/(Et_of_softer_protojet) < overlap_threshold
83 /// the overlapping energy is split between the two protojets;
84 /// otherwise the less energetic protojet is discarded. Called
85 /// ovlim in pxcone.
86 ///
87 /// - pxcone carries out p-scheme recombination, and the resulting
88 /// jets are massless; setting E_scheme_jets = true (default
89 /// false) doesn't change the jet composition, but the final
90 /// momentum sum for the jets is carried out by direct
91 /// four-vector addition instead of p-scheme recombination.
92 ///
93 /// - mode: set to 1 for the e+e- version
94 /// set to 2 for the hadron-hadron version (the default)
95 ///
96 PxConePlugin (double cone_radius_in,
97 double min_jet_energy_in = 5.0,
98 double overlap_threshold_in = 0.5,
99 bool E_scheme_jets_in = false,
100 int mode = 2) :
101 _cone_radius (cone_radius_in ),
102 _min_jet_energy (min_jet_energy_in ),
103 _overlap_threshold (overlap_threshold_in),
104 _E_scheme_jets (E_scheme_jets_in ),
105 _mode (mode ){}
106
107
108 // some functions to return info about parameters ----------------
109
110 /// the cone radius
111 double cone_radius () const {return _cone_radius ;}
112
113 /// minimum jet energy (protojets below this are thrown own before
114 /// merging/splitting) -- called epslon in pxcone
115 double min_jet_energy () const {return _min_jet_energy ;}
116
117 /// Maximum fraction of overlap energy in a jet -- called ovlim in pxcone.
118 double overlap_threshold () const {return _overlap_threshold ;}
119
120 /// if true then the final jets are returned as the E-scheme recombination
121 /// of the particle momenta (by default, pxcone returns massless jets with
122 /// a mean phi,eta type of recombination); regardless of what is
123 /// returned, the internal pxcone jet-finding procedure is
124 /// unaffected.
125 bool E_scheme_jets() const {return _E_scheme_jets ;}
126
127 int mode() const {return _mode ;}
128
129 // the things that are required by base class
130 virtual std::string description () const;
131 virtual void run_clustering(ClusterSequence &) const;
132 /// the plugin mechanism's standard way of accessing the jet radius
133 virtual double R() const {return cone_radius();}
134
135private:
136
137 double _cone_radius ;
138 double _min_jet_energy ;
139 double _overlap_threshold ;
140
141 bool _E_scheme_jets;
142
143 static thread_safety_helpers::FirstTimeTrue _first_time;
144 int _mode; // 1 = e+e-, 2 = hh (default)
145
146 /// print a banner for reference to the 3rd-party code
147 void _print_banner(std::ostream *ostr) const;
148};
149
150FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
151
152#endif // __PXCONEPLUGIN_HH__
a class that allows a user to introduce their own "plugin" jet finder
Implementation of the PxCone algorithm (plugin for fastjet v2.1 upwards)
Definition: PxConePlugin.hh:69
double overlap_threshold() const
Maximum fraction of overlap energy in a jet – called ovlim in pxcone.
bool E_scheme_jets() const
if true then the final jets are returned as the E-scheme recombination of the particle momenta (by de...
PxConePlugin(double cone_radius_in, double min_jet_energy_in=5.0, double overlap_threshold_in=0.5, bool E_scheme_jets_in=false, int mode=2)
constructor for the PxConePlugin, whose arguments have the following meaning:
Definition: PxConePlugin.hh:96
double cone_radius() const
the cone radius
double min_jet_energy() const
minimum jet energy (protojets below this are thrown own before merging/splitting) – called epslon in ...
virtual double R() const
the plugin mechanism's standard way of accessing the jet radius
provides an object wich will return "true" the first time () is called and false afterwards