ATLAS Offline Software
JetPruner.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // JetPruner.cxx
6 
7 #include "JetRec/JetPruner.h"
8 #include <iomanip>
9 #include "fastjet/PseudoJet.hh"
10 #include "fastjet/tools/Pruner.hh"
11 #include "JetEDM/PseudoJetVector.h"
12 #include "JetEDM/FastJetUtils.h"
13 #include "JetRec/JetDumper.h"
14 
15 using std::setw;
16 using fastjet::PseudoJet;
17 using xAOD::JetContainer;
18 
19 //**********************************************************************
20 
21 JetPruner::JetPruner(const std::string& name)
22  : AsgTool(name), m_bld("",this) {
23  declareProperty("JetAlgorithm", m_jetalg ="CamKt");
24  declareProperty("RCut", m_rcut =0.3);
25  declareProperty("ZCut", m_zcut =0.10);
26  declareProperty("JetBuilder", m_bld);
27 }
28 
29 //**********************************************************************
30 
31 JetPruner::~JetPruner() = default;
32 
33 //**********************************************************************
34 
39  ATH_MSG_ERROR("Invalid jet algorithm name: " << m_jetalg);
40  ATH_MSG_ERROR("Allowed values are Kt, CamKt, AntiKt, etc.");
41  return StatusCode::FAILURE;
42  }
43  if ( m_rcut < 0.0 || m_rcut > 10.0 ) {
44  ATH_MSG_ERROR("Invalid value for RCut " << m_rcut);
45  return StatusCode::FAILURE;
46  }
47  if ( m_zcut < 0.0 || m_zcut > 1.0 ) {
48  ATH_MSG_ERROR("Invalid value for ZCut " << m_zcut);
49  return StatusCode::FAILURE;
50  }
51  if ( m_bld.empty() ) {
52  ATH_MSG_ERROR("Unable to retrieve jet builder.");
53  return StatusCode::FAILURE;
54  }
55  return StatusCode::SUCCESS;
56 }
57 
58 //**********************************************************************
59 
60 int JetPruner::groom(const xAOD::Jet& jin,
61  const PseudoJetContainer& pjContainer,
62  xAOD::JetContainer& jets) const {
63  if ( pseudojetRetriever() == nullptr ) {
64  ATH_MSG_WARNING("Pseudojet retriever is null.");
65  return 1;
66  }
67  const PseudoJet* ppjin = pseudojetRetriever()->pseudojet(jin);
68  if ( ppjin == nullptr ) {
69  ATH_MSG_WARNING("Jet does not have a pseudojet.");
70  return 1;
71  }
72  // Prune.
73  fastjet::Pruner pruner(m_fjalg, m_zcut, m_rcut);
74  PseudoJet pjprun = pruner(*ppjin);
75  ATH_MSG_VERBOSE(" Input cluster sequence: " << ppjin->associated_cluster_sequence());
76  ATH_MSG_VERBOSE(" Pruned cluster sequence: " << pjprun.associated_cluster_sequence());
77  // Add jet to collection.
78  xAOD::Jet* pjet = m_bld->add(pjprun, pjContainer, jets, &jin);
79  pjet->setAttribute<float>("RCut", m_rcut);
80  pjet->setAttribute<float>("ZCut", m_zcut);
81  pjet->setAttribute<int>("TransformType", xAOD::JetTransform::Prune);
82  ATH_MSG_DEBUG("Properties after pruning:");
83  ATH_MSG_DEBUG(" ncon: " << pjprun.constituents().size() << "/"
84  << ppjin->constituents().size());
85  ATH_MSG_DEBUG(" nsub: " << pjprun.pieces().size());
86  return 0;
87 }
88 
89 //**********************************************************************
90 
91 void JetPruner::print() const {
92  ATH_MSG_INFO(" Jet algorithm: " << m_jetalg);
93  ATH_MSG_INFO(" R cut factor: " << m_rcut);
94  ATH_MSG_INFO(" Z cut: " << m_zcut);
95  ATH_MSG_INFO(" Jet builder: " << m_bld.name());
96 }
97 
98 //**********************************************************************
JetPruner::m_rcut
float m_rcut
Definition: JetPruner.h:50
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
JetPruner::m_fjalg
fastjet::JetAlgorithm m_fjalg
Definition: JetPruner.h:54
AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
PseudoJetVector.h
IJetGroomer::pseudojetRetriever
virtual const IJetPseudojetRetriever * pseudojetRetriever() const
Return the pseudojet retriever associated with this tool.
Definition: IJetGroomer.cxx:21
JetDumper.h
PseudoJetContainer
Definition: PseudoJetContainer.h:48
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
JetPruner::m_jetalg
std::string m_jetalg
Definition: JetPruner.h:48
JetPruner::~JetPruner
~JetPruner()
JetPruner::print
void print() const
Print the state of the tool.
Definition: JetPruner.cxx:91
JetPruner::m_bld
ToolHandle< IJetFromPseudojet > m_bld
Definition: JetPruner.h:51
JetPruner::initialize
StatusCode initialize()
Dummy implementation of the initialisation function.
Definition: JetPruner.cxx:35
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
xAOD::JetAlgorithmType::fastJetDef
fastjet::JetAlgorithm fastJetDef(ID id)
Definition: FastJetUtils.cxx:20
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
JetPruner::JetPruner
JetPruner(const std::string &name)
Definition: JetPruner.cxx:21
JetPruner.h
xAOD::JetAlgorithmType::ID
ID
//////////////////////////////////////// JetAlgorithmType::ID defines most common physics jet finding...
Definition: JetContainerInfo.h:29
xAOD::JetAlgorithmType::undefined_jet_algorithm
@ undefined_jet_algorithm
Definition: JetContainerInfo.h:40
DataVector
Derived DataVector<T>.
Definition: DataVector.h:794
xAOD::Jet_v1::setAttribute
void setAttribute(const std::string &name, const T &v)
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
IJetPseudojetRetriever::pseudojet
virtual const fastjet::PseudoJet * pseudojet(const xAOD::Jet &jet) const =0
Retrieve the pseudojet associate with a jet.
JetPruner::groom
int groom(const xAOD::Jet &jin, const PseudoJetContainer &, xAOD::JetContainer &jout) const
Transform jet.
Definition: JetPruner.cxx:60
xAOD::Jet_v1
Class describing a jet.
Definition: Jet_v1.h:57
JetPruner::m_zcut
float m_zcut
Definition: JetPruner.h:49
xAOD::JetTransform::Prune
@ Prune
Definition: JetContainerInfo.h:116
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
xAOD::JetAlgorithmType::algId
ID algId(const std::string &n)
Converts a string into a JetAlgorithmType::ID.
Definition: JetContainerInfo.cxx:76
defineDB.jets
list jets
Definition: JetTagCalibration/share/defineDB.py:24
xAOD::JetContainer
JetContainer_v1 JetContainer
Definition of the current "jet container version".
Definition: JetContainer.h:17
FastJetUtils.h