ATLAS Offline Software
JetFinder.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // JetFinder.h
6 
7 #ifndef JetFinder_H
8 #define JetFinder_H
9 
10 
11 // Tool to find jets using fastjet.
12 // The output collection is of type xAOD::JetContainer.
13 
14 #include "AsgTools/AsgTool.h"
17 #include "fastjet/PseudoJet.hh"
18 #include "fastjet/JetDefinition.hh"
19 #include "xAODJet/JetContainer.h"
20 #include "AsgTools/ToolHandle.h"
25 
26 namespace fastjet {
27  class ClusterSequence;
28 }
29 
44 
45 
46 class JetFinder
47 : public asg::AsgTool,
48  virtual public IJetFinder {
50 
51 public:
52 
53  JetFinder(const std::string& name);
54 
55  virtual StatusCode initialize() override;
56 
57  // Find jets and put them in a container. Save ClusterSequence to evtstore
58  virtual int find(const PseudoJetContainer& cont,
59  xAOD::JetContainer & finalJets,
60  xAOD::JetInput::Type contype ) const override;
61 
62  // Find jets and put them in a container - no writes to evet store
63  virtual int findNoSave(const PseudoJetContainer& cont,
64  xAOD::JetContainer & finalJets,
65  xAOD::JetInput::Type contype,
66  fastjet::ClusterSequence*& cs) const override;
67 
68  // Save a cluster sequence in the event store.
69  // So it lives as long as this event.
70  void save(fastjet::ClusterSequence* pcs) const;
71 
72  // Is variable-finding done?
73  bool isVariableR() const;
74 
75  // Dump to log.
76  virtual void print() const override;
77 
78 private: //data
82 
83  // Job options.
84  std::string m_jetalg; // Kt, AntiKt, CamKt
85  float m_jetrad; // Jet size parameter.
86  float m_minrad; // Variable-R min radius
87  float m_massscale; // Variable-R mass scale
88  float m_ptmin; // pT min in MeV
89  float m_ghostarea; // Area for ghosts. 0==>no ghosts.
90  int m_ranopt; // Rand option: 0=fj default, 1=run/event
91 
92  ToolHandle<IJetFromPseudojet> m_bld;
93 
94  // Data
95  fastjet::JetAlgorithm m_fjalg;
97 
98  // Implmentation of find(), findNoSave()
99  int _find(const PseudoJetContainer& cont,
100  xAOD::JetContainer & finalJets,
101  xAOD::JetInput::Type contype,
102  bool doSave,
103  fastjet::ClusterSequence*&) const;
104 
105 };
106 
107 #endif
JetFinder
Tool to find jets.
Definition: JetFinder.h:48
asg::AsgTool
Base class for the dual-use tool implementation classes.
Definition: AsgTool.h:47
fastjet
Definition: FastJetLinkBase.h:22
JetFinder::JetFinder
JetFinder(const std::string &name)
Definition: JetFinder.cxx:37
IJetFinder
Definition: IJetFinder.h:30
JetFinder::m_fjalg
fastjet::JetAlgorithm m_fjalg
Definition: JetFinder.h:95
JetFinder::m_cnameRKey
SG::ReadHandleKey< fastjet::ClusterSequence > m_cnameRKey
Definition: JetFinder.h:80
PseudoJetContainer
Definition: PseudoJetContainer.h:48
IJetFinder.h
SG::ReadHandleKey< xAOD::EventInfo >
IJetFromPseudojet.h
JetFinder::isVariableR
bool isVariableR() const
Definition: JetFinder.cxx:264
JetFinder::find
virtual int find(const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype) const override
Method to find jets from a vector of pseudojet inputs.
Definition: JetFinder.cxx:98
SG::WriteHandleKey< fastjet::ClusterSequence >
jet::ClusterSequence
fastjet::ClusterSequence ClusterSequence
Definition: ClusterSequence.h:21
JetFinder::m_jetrad
float m_jetrad
Definition: JetFinder.h:85
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ReadHandleKey.h
Property holding a SG store/key/clid from which a ReadHandle is made.
JetFinder::m_jetalg
std::string m_jetalg
Definition: JetFinder.h:84
JetFinder::m_minrad
float m_minrad
Definition: JetFinder.h:86
WriteHandleKey.h
Property holding a SG store/key/clid from which a WriteHandle is made.
JetFinder::m_ghostarea
float m_ghostarea
Definition: JetFinder.h:89
xAOD::JetInput::Type
Type
Definition: JetContainerInfo.h:54
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
JetFinder::m_eventinfokey
SG::ReadHandleKey< xAOD::EventInfo > m_eventinfokey
Definition: JetFinder.h:79
JetFinder::initialize
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Definition: JetFinder.cxx:54
JetFinder::m_massscale
float m_massscale
Definition: JetFinder.h:87
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
JetFinder::m_ptmin
float m_ptmin
Definition: JetFinder.h:88
EventInfo.h
JetFinder::m_ranopt
int m_ranopt
Definition: JetFinder.h:90
JetFinder::save
void save(fastjet::ClusterSequence *pcs) const
Definition: JetFinder.cxx:242
PseudoJetContainer.h
JetContainer.h
ASG_TOOL_CLASS
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
Definition: AsgToolMacros.h:68
JetFinder::_find
int _find(const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, bool doSave, fastjet::ClusterSequence *&) const
Definition: JetFinder.cxx:119
JetFinder::findNoSave
virtual int findNoSave(const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, fastjet::ClusterSequence *&cs) const override
Definition: JetFinder.cxx:109
ToolHandle.h
AsgTool.h
python.CreateTierZeroArgdict.pcs
pcs
Definition: CreateTierZeroArgdict.py:200
JetFinder::print
virtual void print() const override
Print the state of the tool.
Definition: JetFinder.cxx:270
JetFinder::m_cnameWKey
SG::WriteHandleKey< fastjet::ClusterSequence > m_cnameWKey
Definition: JetFinder.h:81
JetFinder::m_isVariableR
bool m_isVariableR
Definition: JetFinder.h:96
JetFinder::m_bld
ToolHandle< IJetFromPseudojet > m_bld
Definition: JetFinder.h:92