ATLAS Offline Software
Loading...
Searching...
No Matches
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"
20#include "AsgTools/ToolHandle.h"
25
26namespace fastjet {
27 class ClusterSequence;
28}
29
44
45
47: public asg::AsgTool,
48 virtual public IJetFinder {
50
51public:
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,
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
78private: //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
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
Property holding a SG store/key/clid from which a ReadHandle is made.
Property holding a SG store/key/clid from which a WriteHandle is made.
IJetFinder is a dual-use tool interface for for a tool that modifies a jet collection.
Definition IJetFinder.h:30
std::string m_jetalg
Definition JetFinder.h:84
float m_massscale
Definition JetFinder.h:87
int m_ranopt
Definition JetFinder.h:90
virtual int findNoSave(const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, fastjet::ClusterSequence *&cs) const override
virtual void print() const override
Print the state of the tool.
int _find(const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, bool doSave, fastjet::ClusterSequence *&) const
SG::WriteHandleKey< fastjet::ClusterSequence > m_cnameWKey
Definition JetFinder.h:81
float m_ptmin
Definition JetFinder.h:88
float m_ghostarea
Definition JetFinder.h:89
float m_jetrad
Definition JetFinder.h:85
bool m_isVariableR
Definition JetFinder.h:96
bool isVariableR() const
JetFinder(const std::string &name)
Definition JetFinder.cxx:37
void save(fastjet::ClusterSequence *pcs) const
SG::ReadHandleKey< fastjet::ClusterSequence > m_cnameRKey
Definition JetFinder.h:80
float m_minrad
Definition JetFinder.h:86
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Definition JetFinder.cxx:54
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
ToolHandle< IJetFromPseudojet > m_bld
Definition JetFinder.h:92
SG::ReadHandleKey< xAOD::EventInfo > m_eventinfokey
Definition JetFinder.h:79
fastjet::JetAlgorithm m_fjalg
Definition JetFinder.h:95
Property holding a SG store/key/clid from which a ReadHandle is made.
Property holding a SG store/key/clid from which a WriteHandle is made.
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
fastjet::ClusterSequence ClusterSequence
JetContainer_v1 JetContainer
Definition of the current "jet container version".