Fork me on GitHub

source: git/external/fastjet/ClusterSequenceStructure.hh@ 41376ac

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 41376ac was d69dfe4, checked in by pavel <pavel@…>, 11 years ago

upgrade fastjet to 3.0.6

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