ATLAS Offline Software
Loading...
Searching...
No Matches
JetClusterer.h
Go to the documentation of this file.
1// this file is -*- C++ -*-
2/*
3 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
4*/
5
6#ifndef JETREC_JETCLUSTERER_H
7#define JETREC_JETCLUSTERER_H
18
19
24
26#include "AsgTools/AsgTool.h"
28
34
35#include "fastjet/PseudoJet.hh"
36#include "fastjet/AreaDefinition.hh"
37#include "fastjet/JetDefinition.hh"
38
41
42
44: public asg::AsgTool,
45 virtual public JetProvider<xAOD::JetAuxContainer> {
47
48public:
49
50 // Can't use "using ctor" because of incompatiblity with pyroot in AnalysisBase
51 JetClusterer(const std::string &name)
52 : AsgTool(name),
53 m_jetRankAccessor("jetRank")
54 {}
55
56 StatusCode initialize() override;
57
60 std::pair<std::unique_ptr<xAOD::JetContainer>, std::unique_ptr<SG::IAuxStore> > getJets() const override;
61
62
63protected:
64
66 fastjet::AreaDefinition buildAreaDefinition(bool &seedsok) const ;
67
69 std::unique_ptr<fastjet::ClusterSequence> buildClusterSequence(const PseudoJetVector* pseudoJetvector) const ;
70
72 void processPseudoJet(const fastjet::PseudoJet& pj, const PseudoJetContainer& pjCont, xAOD::JetContainer* jets, const xAOD::Vertex* vertex) const;
73
76
78 SG::ReadHandleKey<PseudoJetContainer> m_inputPseudoJets {this, "InputPseudoJets", "inputpseudojet", "input constituents"};
79
81 SG::WriteHandleKey<PseudoJetVector> m_finalPseudoJets {this, "FinalPseudoJets_DONOTSET", "", "output pseudojets -- autoconfigured name"};
82 SG::WriteHandleKey<jet::ClusterSequence> m_clusterSequence {this, "ClusterSequence_DONOTSET", "", "output pseudojets -- autoconfigured name"};
83 Gaudi::Property<std::string> m_jetRank {this, "jetRankName", "jetRank", "name for accessor for jet rank"};
84
85 // Job options.
86 Gaudi::Property<std::string> m_jetalg {this, "JetAlgorithm", "AntiKt", "alg type : AntiKt, Kt, CA..."};
87 Gaudi::Property<float> m_jetrad {this, "JetRadius", 0.4 , "jet size parameter"};
88 Gaudi::Property<float> m_ptmin {this, "PtMin", 0.0, "pT min in MeV"};
89 Gaudi::Property<float> m_ghostarea {this, "GhostArea", 0.0, "Area for ghosts. 0==>no ghosts."};
90 Gaudi::Property<int> m_ranopt {this, "RandomOption", 1, "Rand option: 0=fj default, 1=run/event"};
91
92 Gaudi::Property<int> m_inputType {this, "JetInputType", 0, "input type as in xAOD::JetInput (see xAODJet/JetContainerInfo.h)"};
93
94
95 // internal data set from properties during initialize()
96 fastjet::JetAlgorithm m_fjalg;
97 bool m_useArea{};
98
99 Gaudi::Property<float> m_minrad {this, "VariableRMinRadius", -1.0, "Variable-R min radius" };
100 Gaudi::Property<float> m_massscale {this, "VariableRMassScale", -1.0, "Variable-R mass scale" };
102 bool isVariableR() const { return m_isVariableR;}
103
105
106};
107
108
110
112 void seedsFromEventInfo(const xAOD::EventInfo* ei, std::vector<int> &seeds);
113
114}
115
116#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.
Base class for elements of a container that can have aux data.
std::vector< fastjet::PseudoJet > PseudoJetVector
Gaudi::Property< float > m_ghostarea
bool isVariableR() const
fastjet::AreaDefinition buildAreaDefinition(bool &seedsok) const
Build the area definition when running with area calculation. The seedsok flag is set to false when e...
std::pair< std::unique_ptr< xAOD::JetContainer >, std::unique_ptr< SG::IAuxStore > > getJets() const override
Return the final jets with their aux store.
Gaudi::Property< std::string > m_jetalg
SG::WriteHandleKey< jet::ClusterSequence > m_clusterSequence
Gaudi::Property< int > m_inputType
Gaudi::Property< float > m_massscale
Gaudi::Property< float > m_ptmin
Gaudi::Property< std::string > m_jetRank
SG::ReadHandleKey< xAOD::EventInfo > m_eventinfokey
Handle to EventInfo. This is used to get the evt&run number to set the Ghost area random seeds.
Gaudi::Property< int > m_ranopt
StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< float > m_minrad
void processPseudoJet(const fastjet::PseudoJet &pj, const PseudoJetContainer &pjCont, xAOD::JetContainer *jets, const xAOD::Vertex *vertex) const
translate to xAOD::Jet
SG::AuxElement::Accessor< int > m_jetRankAccessor
Gaudi::Property< float > m_jetrad
std::unique_ptr< fastjet::ClusterSequence > buildClusterSequence(const PseudoJetVector *pseudoJetvector) const
Used to create the cluster sequence.
fastjet::JetAlgorithm m_fjalg
SG::ReadHandleKey< PseudoJetContainer > m_inputPseudoJets
Handle Input PseudoJetContainer.
SG::WriteHandleKey< PseudoJetVector > m_finalPseudoJets
used to build the key under which the final PJ will be stored in evtStore()
JetClusterer(const std::string &name)
Concrete class that implements the recording of jets & aux container to StoreGate via an externally p...
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
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
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
void seedsFromEventInfo(const xAOD::EventInfo *ei, std::vector< int > &seeds)
Fill seeds vector from run & event number. This functio is separated from the class so it's easier to...
EventInfo_v1 EventInfo
Definition of the latest event info version.
Vertex_v1 Vertex
Define the latest version of the vertex class.
JetContainer_v1 JetContainer
Definition of the current "jet container version".