|
||||
Warning, file /include/fastjet/ClusterSequenceStructure.hh was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001 //FJSTARTHEADER 0002 // $Id$ 0003 // 0004 // Copyright (c) 2005-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez 0005 // 0006 //---------------------------------------------------------------------- 0007 // This file is part of FastJet. 0008 // 0009 // FastJet is free software; you can redistribute it and/or modify 0010 // it under the terms of the GNU General Public License as published by 0011 // the Free Software Foundation; either version 2 of the License, or 0012 // (at your option) any later version. 0013 // 0014 // The algorithms that underlie FastJet have required considerable 0015 // development. They are described in the original FastJet paper, 0016 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use 0017 // FastJet as part of work towards a scientific publication, please 0018 // quote the version you use and include a citation to the manual and 0019 // optionally also to hep-ph/0512210. 0020 // 0021 // FastJet is distributed in the hope that it will be useful, 0022 // but WITHOUT ANY WARRANTY; without even the implied warranty of 0023 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 0024 // GNU General Public License for more details. 0025 // 0026 // You should have received a copy of the GNU General Public License 0027 // along with FastJet. If not, see <http://www.gnu.org/licenses/>. 0028 //---------------------------------------------------------------------- 0029 //FJENDHEADER 0030 0031 0032 #ifndef __FASTJET_CLUSTER_SEQUENCE_STRUCTURE_HH__ 0033 #define __FASTJET_CLUSTER_SEQUENCE_STRUCTURE_HH__ 0034 0035 #include "fastjet/internal/base.hh" 0036 #include "fastjet/SharedPtr.hh" 0037 #include "fastjet/PseudoJetStructureBase.hh" 0038 0039 #include <vector> 0040 0041 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 0042 0043 /// @ingroup extra_info 0044 /// \class ClusterSequenceStructure 0045 /// 0046 /// Contains any information related to the clustering that should be 0047 /// directly accessible to PseudoJet. 0048 /// 0049 /// By default, this class implements basic access to the 0050 /// ClusterSequence related to a PseudoJet (like its constituents or 0051 /// its area). But it can be overloaded in order e.g. to give access 0052 /// to the jet substructure. 0053 /// 0054 // Design question: Do we only put the methods that can be overloaded 0055 // or do we put everything a PJ can have access to? I think both cost 0056 // the same number of indirections. The first option limits the amount 0057 // of coding and maybe has a clearer structure. The second is more 0058 // consistent (everything related to the same thing is at the same 0059 // place) and gives better access for derived classes. We'll go for 0060 // the second option. 0061 class ClusterSequenceStructure : public PseudoJetStructureBase{ 0062 public: 0063 /// default ctor 0064 ClusterSequenceStructure() : _associated_cs(NULL){} 0065 0066 /// ctor with initialisation to a given ClusterSequence 0067 /// 0068 /// In principle, this is reserved for initialisation by the parent 0069 /// ClusterSequence 0070 ClusterSequenceStructure(const ClusterSequence *cs){ 0071 set_associated_cs(cs); 0072 }; 0073 0074 /// default (virtual) dtor 0075 virtual ~ClusterSequenceStructure(); 0076 0077 /// description 0078 virtual std::string description() const FASTJET_OVERRIDE{ 0079 return "PseudoJet with an associated ClusterSequence"; 0080 } 0081 0082 //------------------------------------------------------------- 0083 /// @name Direct access to the associated ClusterSequence object. 0084 /// 0085 /// Get access to the associated ClusterSequence (if any) 0086 //\{ 0087 //------------------------------------------------------------- 0088 /// returns true if there is an associated ClusterSequence 0089 virtual bool has_associated_cluster_sequence() const FASTJET_OVERRIDE{ return true;} 0090 0091 /// get a (const) pointer to the parent ClusterSequence (NULL if 0092 /// inexistent) 0093 virtual const ClusterSequence* associated_cluster_sequence() const FASTJET_OVERRIDE; 0094 0095 /// returns true if there is a valid associated ClusterSequence 0096 virtual bool has_valid_cluster_sequence() const FASTJET_OVERRIDE; 0097 0098 /// if the jet has a valid associated cluster sequence then return a 0099 /// pointer to it; otherwise throw an error 0100 virtual const ClusterSequence * validated_cs() const FASTJET_OVERRIDE; 0101 0102 #ifndef __FJCORE__ 0103 /// if the jet has valid area information then return a pointer to 0104 /// the associated ClusterSequenceAreaBase object; otherwise throw an error 0105 virtual const ClusterSequenceAreaBase * validated_csab() const FASTJET_OVERRIDE; 0106 #endif // __FJCORE__ 0107 0108 /// set the associated csw 0109 virtual void set_associated_cs(const ClusterSequence * new_cs){ 0110 _associated_cs = new_cs; 0111 } 0112 //\} 0113 0114 //------------------------------------------------------------- 0115 /// @name Methods for access to information about jet structure 0116 /// 0117 /// These allow access to jet constituents, and other jet 0118 /// subtructure information. They only work if the jet is associated 0119 /// with a ClusterSequence. 0120 //------------------------------------------------------------- 0121 //\{ 0122 0123 /// check if it has been recombined with another PseudoJet in which 0124 /// case, return its partner through the argument. Otherwise, 0125 /// 'partner' is set to 0. 0126 /// 0127 /// an Error is thrown if this PseudoJet has no currently valid 0128 /// associated ClusterSequence 0129 virtual bool has_partner(const PseudoJet &reference, PseudoJet &partner) const FASTJET_OVERRIDE; 0130 0131 /// check if it has been recombined with another PseudoJet in which 0132 /// case, return its child through the argument. Otherwise, 'child' 0133 /// is set to 0. 0134 /// 0135 /// an Error is thrown if this PseudoJet has no currently valid 0136 /// associated ClusterSequence 0137 virtual bool has_child(const PseudoJet &reference, PseudoJet &child) const FASTJET_OVERRIDE; 0138 0139 /// check if it is the product of a recombination, in which case 0140 /// return the 2 parents through the 'parent1' and 'parent2' 0141 /// arguments. Otherwise, set these to 0. 0142 /// 0143 /// an Error is thrown if this PseudoJet has no currently valid 0144 /// associated ClusterSequence 0145 virtual bool has_parents(const PseudoJet &reference, PseudoJet &parent1, PseudoJet &parent2) const FASTJET_OVERRIDE; 0146 0147 /// check if the reference PseudoJet is contained in the second one 0148 /// passed as argument. 0149 /// 0150 /// an Error is thrown if this PseudoJet has no currently valid 0151 /// associated ClusterSequence 0152 /// 0153 /// false is returned if the 2 PseudoJet do not belong the same 0154 /// ClusterSequence 0155 virtual bool object_in_jet(const PseudoJet &reference, const PseudoJet &jet) const FASTJET_OVERRIDE; 0156 0157 /// return true if the structure supports constituents. 0158 /// 0159 /// an Error is thrown if this PseudoJet has no currently valid 0160 /// associated ClusterSequence 0161 virtual bool has_constituents() const FASTJET_OVERRIDE; 0162 0163 /// retrieve the constituents. 0164 /// 0165 /// an Error is thrown if this PseudoJet has no currently valid 0166 /// associated ClusterSequence 0167 virtual std::vector<PseudoJet> constituents(const PseudoJet &reference) const FASTJET_OVERRIDE; 0168 0169 0170 /// return true if the structure supports exclusive_subjets. 0171 /// 0172 /// an Error is thrown if this PseudoJet has no currently valid 0173 /// associated ClusterSequence 0174 virtual bool has_exclusive_subjets() const FASTJET_OVERRIDE; 0175 0176 /// return a vector of all subjets of the current jet (in the sense 0177 /// of the exclusive algorithm) that would be obtained when running 0178 /// the algorithm with the given dcut. 0179 /// 0180 /// Time taken is O(m ln m), where m is the number of subjets that 0181 /// are found. If m gets to be of order of the total number of 0182 /// constituents in the jet, this could be substantially slower than 0183 /// just getting that list of constituents. 0184 /// 0185 /// an Error is thrown if this PseudoJet has no currently valid 0186 /// associated ClusterSequence 0187 virtual std::vector<PseudoJet> exclusive_subjets(const PseudoJet &reference, const double & dcut) const FASTJET_OVERRIDE; 0188 0189 /// return the size of exclusive_subjets(...); still n ln n with same 0190 /// coefficient, but marginally more efficient than manually taking 0191 /// exclusive_subjets.size() 0192 /// 0193 /// an Error is thrown if this PseudoJet has no currently valid 0194 /// associated ClusterSequence 0195 virtual int n_exclusive_subjets(const PseudoJet &reference, const double & dcut) const FASTJET_OVERRIDE; 0196 0197 /// return the list of subjets obtained by unclustering the supplied 0198 /// jet down to nsub subjets (or all constituents if there are fewer 0199 /// than nsub). 0200 /// 0201 /// requires nsub ln nsub time 0202 /// 0203 /// an Error is thrown if this PseudoJet has no currently valid 0204 /// associated ClusterSequence 0205 virtual std::vector<PseudoJet> exclusive_subjets_up_to (const PseudoJet &reference, int nsub) const FASTJET_OVERRIDE; 0206 0207 /// return the dij that was present in the merging nsub+1 -> nsub 0208 /// subjets inside this jet. 0209 /// 0210 /// an Error is thrown if this PseudoJet has no currently valid 0211 /// associated ClusterSequence 0212 virtual double exclusive_subdmerge(const PseudoJet &reference, int nsub) const FASTJET_OVERRIDE; 0213 0214 /// return the maximum dij that occurred in the whole event at the 0215 /// stage that the nsub+1 -> nsub merge of subjets occurred inside 0216 /// this jet. 0217 /// 0218 /// an Error is thrown if this PseudoJet has no currently valid 0219 /// associated ClusterSequence 0220 virtual double exclusive_subdmerge_max(const PseudoJet &reference, int nsub) const FASTJET_OVERRIDE; 0221 0222 0223 //------------------------------------------------------------------- 0224 // information related to the pieces of the jet 0225 //------------------------------------------------------------------- 0226 /// by convention, a jet associated with a ClusterSequence will have 0227 /// its parents as pieces 0228 virtual bool has_pieces(const PseudoJet &reference) const FASTJET_OVERRIDE; 0229 0230 /// by convention, a jet associated with a ClusterSequence will have 0231 /// its parents as pieces 0232 /// 0233 /// if it has no parents, then there will only be a single piece: 0234 /// itself 0235 /// 0236 /// Note that to answer that question, we need to access the cluster 0237 /// sequence. If the cluster sequence has gone out of scope, an 0238 /// error will be thrown 0239 virtual std::vector<PseudoJet> pieces(const PseudoJet &reference) const FASTJET_OVERRIDE; 0240 0241 0242 // the following ones require a computation of the area in the 0243 // parent ClusterSequence (See ClusterSequenceAreaBase for details) 0244 //------------------------------------------------------------------ 0245 #ifndef __FJCORE__ 0246 0247 /// check if it has a defined area 0248 virtual bool has_area() const FASTJET_OVERRIDE; 0249 0250 /// return the jet (scalar) area. 0251 /// throws an Error if there is no support for area in the parent CS 0252 virtual double area(const PseudoJet &reference) const FASTJET_OVERRIDE; 0253 0254 /// return the error (uncertainty) associated with the determination 0255 /// of the area of this jet. 0256 /// throws an Error if there is no support for area in the parent CS 0257 virtual double area_error(const PseudoJet &reference) const FASTJET_OVERRIDE; 0258 0259 /// return the jet 4-vector area. 0260 /// throws an Error if there is no support for area in the parent CS 0261 virtual PseudoJet area_4vector(const PseudoJet &reference) const FASTJET_OVERRIDE; 0262 0263 /// true if this jet is made exclusively of ghosts. 0264 /// throws an Error if there is no support for area in the parent CS 0265 virtual bool is_pure_ghost(const PseudoJet &reference) const FASTJET_OVERRIDE; 0266 0267 #endif // __FJCORE__ 0268 //\} --- end of jet structure ------------------------------------- 0269 0270 protected: 0271 const ClusterSequence *_associated_cs; 0272 }; 0273 0274 FASTJET_END_NAMESPACE 0275 0276 #endif // __FASTJET_CLUSTER_SEQUENCE_STRUCTURE_HH__
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |