Fork me on GitHub

source: git/external/fastjet/ClusterSequenceStructure.hh@ 5eb063e

Last change on this file since 5eb063e was b7b836a, checked in by Pavel Demin <pavel-demin@…>, 7 years ago

update FastJet library to 3.3.1 and FastJet Contrib library to 1.036

  • Property mode set to 100644
File size: 11.1 KB
RevLine 
[35cdc46]1//FJSTARTHEADER
[b7b836a]2// $Id: ClusterSequenceStructure.hh 4354 2018-04-22 07:12:37Z salam $
[d7d2da3]3//
[b7b836a]4// Copyright (c) 2005-2018, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
[d7d2da3]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
[35cdc46]15// development. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
[d7d2da3]17// FastJet as part of work towards a scientific publication, please
[35cdc46]18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
[d7d2da3]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//----------------------------------------------------------------------
[35cdc46]29//FJENDHEADER
[d7d2da3]30
31
32#ifndef __FASTJET_CLUSTER_SEQUENCE_STRUCTURE_HH__
33#define __FASTJET_CLUSTER_SEQUENCE_STRUCTURE_HH__
34
35#include "fastjet/internal/base.hh"
36#include "fastjet/SharedPtr.hh"
37#include "fastjet/PseudoJetStructureBase.hh"
38
39#include <vector>
40
41FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
42
43/// @ingroup extra_info
44/// \class ClusterSequenceStructure
45///
46/// Contains any information related to the clustering that should be
47/// directly accessible to PseudoJet.
48///
49/// By default, this class implements basic access to the
50/// ClusterSequence related to a PseudoJet (like its constituents or
51/// its area). But it can be overloaded in order e.g. to give access
52/// to the jet substructure.
53///
54// Design question: Do we only put the methods that can be overloaded
55// or do we put everything a PJ can have access to? I think both cost
56// the same number of indirections. The first option limits the amount
57// of coding and maybe has a clearer structure. The second is more
58// consistent (everything related to the same thing is at the same
59// place) and gives better access for derived classes. We'll go for
60// the second option.
61class ClusterSequenceStructure : public PseudoJetStructureBase{
62public:
63 /// default ctor
64 ClusterSequenceStructure() : _associated_cs(NULL){}
65
66 /// ctor with initialisation to a given ClusterSequence
67 ///
68 /// In principle, this is reserved for initialisation by the parent
69 /// ClusterSequence
70 ClusterSequenceStructure(const ClusterSequence *cs){
71 set_associated_cs(cs);
72 };
73
74 /// default (virtual) dtor
75 virtual ~ClusterSequenceStructure();
76
77 /// description
[1d208a2]78 virtual std::string description() const FASTJET_OVERRIDE{
79 return "PseudoJet with an associated ClusterSequence";
80 }
[d7d2da3]81
82 //-------------------------------------------------------------
83 /// @name Direct access to the associated ClusterSequence object.
84 ///
85 /// Get access to the associated ClusterSequence (if any)
86 //\{
87 //-------------------------------------------------------------
88 /// returns true if there is an associated ClusterSequence
[1d208a2]89 virtual bool has_associated_cluster_sequence() const FASTJET_OVERRIDE{ return true;}
[d7d2da3]90
91 /// get a (const) pointer to the parent ClusterSequence (NULL if
92 /// inexistent)
[1d208a2]93 virtual const ClusterSequence* associated_cluster_sequence() const FASTJET_OVERRIDE;
[d7d2da3]94
95 /// returns true if there is a valid associated ClusterSequence
[1d208a2]96 virtual bool has_valid_cluster_sequence() const FASTJET_OVERRIDE;
[d7d2da3]97
98 /// if the jet has a valid associated cluster sequence then return a
99 /// pointer to it; otherwise throw an error
[1d208a2]100 virtual const ClusterSequence * validated_cs() const FASTJET_OVERRIDE;
[d7d2da3]101
[d69dfe4]102#ifndef __FJCORE__
[d7d2da3]103 /// if the jet has valid area information then return a pointer to
104 /// the associated ClusterSequenceAreaBase object; otherwise throw an error
[1d208a2]105 virtual const ClusterSequenceAreaBase * validated_csab() const FASTJET_OVERRIDE;
[d69dfe4]106#endif // __FJCORE__
[d7d2da3]107
108 /// set the associated csw
109 virtual void set_associated_cs(const ClusterSequence * new_cs){
110 _associated_cs = new_cs;
111 }
112 //\}
113
114 //-------------------------------------------------------------
115 /// @name Methods for access to information about jet structure
116 ///
117 /// These allow access to jet constituents, and other jet
118 /// subtructure information. They only work if the jet is associated
119 /// with a ClusterSequence.
120 //-------------------------------------------------------------
121 //\{
122
123 /// check if it has been recombined with another PseudoJet in which
124 /// case, return its partner through the argument. Otherwise,
125 /// 'partner' is set to 0.
126 ///
127 /// an Error is thrown if this PseudoJet has no currently valid
128 /// associated ClusterSequence
[1d208a2]129 virtual bool has_partner(const PseudoJet &reference, PseudoJet &partner) const FASTJET_OVERRIDE;
[d7d2da3]130
131 /// check if it has been recombined with another PseudoJet in which
132 /// case, return its child through the argument. Otherwise, 'child'
133 /// is set to 0.
134 ///
135 /// an Error is thrown if this PseudoJet has no currently valid
136 /// associated ClusterSequence
[1d208a2]137 virtual bool has_child(const PseudoJet &reference, PseudoJet &child) const FASTJET_OVERRIDE;
[d7d2da3]138
139 /// check if it is the product of a recombination, in which case
140 /// return the 2 parents through the 'parent1' and 'parent2'
141 /// arguments. Otherwise, set these to 0.
142 ///
143 /// an Error is thrown if this PseudoJet has no currently valid
144 /// associated ClusterSequence
[1d208a2]145 virtual bool has_parents(const PseudoJet &reference, PseudoJet &parent1, PseudoJet &parent2) const FASTJET_OVERRIDE;
[d7d2da3]146
147 /// check if the reference PseudoJet is contained in the second one
148 /// passed as argument.
149 ///
150 /// an Error is thrown if this PseudoJet has no currently valid
151 /// associated ClusterSequence
152 ///
153 /// false is returned if the 2 PseudoJet do not belong the same
154 /// ClusterSequence
[1d208a2]155 virtual bool object_in_jet(const PseudoJet &reference, const PseudoJet &jet) const FASTJET_OVERRIDE;
[d7d2da3]156
157 /// return true if the structure supports constituents.
158 ///
159 /// an Error is thrown if this PseudoJet has no currently valid
160 /// associated ClusterSequence
[1d208a2]161 virtual bool has_constituents() const FASTJET_OVERRIDE;
[d7d2da3]162
163 /// retrieve the constituents.
164 ///
165 /// an Error is thrown if this PseudoJet has no currently valid
166 /// associated ClusterSequence
[1d208a2]167 virtual std::vector<PseudoJet> constituents(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]168
169
170 /// return true if the structure supports exclusive_subjets.
171 ///
172 /// an Error is thrown if this PseudoJet has no currently valid
173 /// associated ClusterSequence
[1d208a2]174 virtual bool has_exclusive_subjets() const FASTJET_OVERRIDE;
[d7d2da3]175
176 /// return a vector of all subjets of the current jet (in the sense
177 /// of the exclusive algorithm) that would be obtained when running
178 /// the algorithm with the given dcut.
179 ///
180 /// Time taken is O(m ln m), where m is the number of subjets that
181 /// are found. If m gets to be of order of the total number of
182 /// constituents in the jet, this could be substantially slower than
183 /// just getting that list of constituents.
184 ///
185 /// an Error is thrown if this PseudoJet has no currently valid
186 /// associated ClusterSequence
[1d208a2]187 virtual std::vector<PseudoJet> exclusive_subjets(const PseudoJet &reference, const double & dcut) const FASTJET_OVERRIDE;
[d7d2da3]188
189 /// return the size of exclusive_subjets(...); still n ln n with same
190 /// coefficient, but marginally more efficient than manually taking
191 /// exclusive_subjets.size()
192 ///
193 /// an Error is thrown if this PseudoJet has no currently valid
194 /// associated ClusterSequence
[1d208a2]195 virtual int n_exclusive_subjets(const PseudoJet &reference, const double & dcut) const FASTJET_OVERRIDE;
[d7d2da3]196
197 /// return the list of subjets obtained by unclustering the supplied
198 /// jet down to nsub subjets (or all constituents if there are fewer
199 /// than nsub).
200 ///
201 /// requires nsub ln nsub time
202 ///
203 /// an Error is thrown if this PseudoJet has no currently valid
204 /// associated ClusterSequence
[1d208a2]205 virtual std::vector<PseudoJet> exclusive_subjets_up_to (const PseudoJet &reference, int nsub) const FASTJET_OVERRIDE;
[d7d2da3]206
207 /// return the dij that was present in the merging nsub+1 -> nsub
208 /// subjets inside this jet.
209 ///
210 /// an Error is thrown if this PseudoJet has no currently valid
211 /// associated ClusterSequence
[1d208a2]212 virtual double exclusive_subdmerge(const PseudoJet &reference, int nsub) const FASTJET_OVERRIDE;
[d7d2da3]213
214 /// return the maximum dij that occurred in the whole event at the
215 /// stage that the nsub+1 -> nsub merge of subjets occurred inside
216 /// this jet.
217 ///
218 /// an Error is thrown if this PseudoJet has no currently valid
219 /// associated ClusterSequence
[1d208a2]220 virtual double exclusive_subdmerge_max(const PseudoJet &reference, int nsub) const FASTJET_OVERRIDE;
[d7d2da3]221
222
223 //-------------------------------------------------------------------
224 // information related to the pieces of the jet
225 //-------------------------------------------------------------------
226 /// by convention, a jet associated with a ClusterSequence will have
227 /// its parents as pieces
[1d208a2]228 virtual bool has_pieces(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]229
230 /// by convention, a jet associated with a ClusterSequence will have
231 /// its parents as pieces
232 ///
233 /// if it has no parents, then there will only be a single piece:
234 /// itself
235 ///
236 /// Note that to answer that question, we need to access the cluster
237 /// sequence. If the cluster sequence has gone out of scope, an
238 /// error will be thrown
[1d208a2]239 virtual std::vector<PseudoJet> pieces(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]240
241
242 // the following ones require a computation of the area in the
243 // parent ClusterSequence (See ClusterSequenceAreaBase for details)
244 //------------------------------------------------------------------
[d69dfe4]245#ifndef __FJCORE__
[d7d2da3]246
247 /// check if it has a defined area
[1d208a2]248 virtual bool has_area() const FASTJET_OVERRIDE;
[d7d2da3]249
250 /// return the jet (scalar) area.
251 /// throws an Error if there is no support for area in the parent CS
[1d208a2]252 virtual double area(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]253
254 /// return the error (uncertainty) associated with the determination
255 /// of the area of this jet.
256 /// throws an Error if there is no support for area in the parent CS
[1d208a2]257 virtual double area_error(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]258
259 /// return the jet 4-vector area.
260 /// throws an Error if there is no support for area in the parent CS
[1d208a2]261 virtual PseudoJet area_4vector(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]262
263 /// true if this jet is made exclusively of ghosts.
264 /// throws an Error if there is no support for area in the parent CS
[1d208a2]265 virtual bool is_pure_ghost(const PseudoJet &reference) const FASTJET_OVERRIDE;
[d7d2da3]266
[d69dfe4]267#endif // __FJCORE__
[d7d2da3]268 //\} --- end of jet structure -------------------------------------
269
270protected:
271 const ClusterSequence *_associated_cs;
272};
273
274FASTJET_END_NAMESPACE
275
276#endif // __FASTJET_CLUSTER_SEQUENCE_STRUCTURE_HH__
Note: See TracBrowser for help on using the repository browser.