ATLAS Offline Software
JetTrimmer.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 // JetTrimmer.cxx
6 
7 #include "JetRec/JetTrimmer.h"
8 #include <iomanip>
9 #include "fastjet/PseudoJet.hh"
10 #include "fastjet/JetDefinition.hh"
11 #include "fastjet/Selector.hh"
12 #include "fastjet/tools/Filter.hh"
13 #include "JetEDM/PseudoJetVector.h"
14 
15 using std::setw;
16 using fastjet::PseudoJet;
17 using xAOD::JetContainer;
18 
19 //**********************************************************************
20 
21 JetTrimmer::JetTrimmer(const std::string& name)
22  : AsgTool(name), m_bld("",this) {
23  declareProperty("RClus", m_rclus =0.3);
24  declareProperty("PtFrac", m_ptfrac =0.03);
25  declareProperty("JetBuilder", m_bld);
26 }
27 
28 //**********************************************************************
29 
30 JetTrimmer::~JetTrimmer() = default;
31 
32 //**********************************************************************
33 
35  if ( m_rclus < 0.0 || m_rclus > 10.0 ) {
36  ATH_MSG_ERROR("Invalid value for RClus " << m_rclus);
37  return StatusCode::FAILURE;
38  }
39  if ( m_ptfrac < 0.0 || m_ptfrac > 1.0 ) {
40  ATH_MSG_ERROR("Invalid value for PtFrac " << m_ptfrac);
41  return StatusCode::FAILURE;
42  }
43  if ( m_bld.empty() ) {
44  ATH_MSG_ERROR("Unable to retrieve jet builder.");
45  return StatusCode::FAILURE;
46  }
47  ATH_CHECK( m_bld.retrieve() );
48  return StatusCode::SUCCESS;
49 }
50 
51 //**********************************************************************
52 
54  const PseudoJetContainer& pjContainer,
55  xAOD::JetContainer& jets) const {
56  if ( pseudojetRetriever() == nullptr ) {
57  ATH_MSG_WARNING("Pseudojet retriever is null.");
58  return 1;
59  }
60  const PseudoJet* ppjin = pseudojetRetriever()->pseudojet(jin);
61  if ( ppjin == nullptr ) {
62  ATH_MSG_WARNING("Jet does not have a pseudojet.");
63  return 1;
64  }
65  // Trim.
66  fastjet::Filter trimmer(fastjet::JetDefinition(fastjet::kt_algorithm, m_rclus),
67  fastjet::SelectorPtFractionMin(m_ptfrac));
68  PseudoJet pjtrim = trimmer(*ppjin);
69  ATH_MSG_VERBOSE(" Input cluster sequence: " << ppjin->associated_cluster_sequence());
70  ATH_MSG_VERBOSE(" Trimmed cluster sequence: " << pjtrim.associated_cluster_sequence());
71  int nptrim = pjtrim.pieces().size();
72  // Add jet to collection.
73  xAOD::Jet* pjet = m_bld->add(pjtrim, pjContainer, jets, &jin);
74  if ( pjet == nullptr ) {
75  ATH_MSG_ERROR("Unable to add jet to container");
76  return 1;
77  }
78  pjet->setAttribute<int>("TransformType", xAOD::JetTransform::Trim);
79  pjet->setAttribute("RClus", float(m_rclus));
80  pjet->setAttribute("PtFrac", float(m_ptfrac));
81  pjet->setAttribute<int>("NTrimSubjets", nptrim);
82  ATH_MSG_DEBUG("Properties after trimming:");
83  ATH_MSG_DEBUG(" ncon: " << pjtrim.constituents().size() << "/"
84  << ppjin->constituents().size());
85  ATH_MSG_DEBUG(" nsub: " << nptrim);
86 
87  ATH_MSG_DEBUG("Added jet to container.");
88 
89  return 0;
90 }
91 
92 //**********************************************************************
93 
94 void JetTrimmer::print() const {
95  ATH_MSG_INFO(" Recluster R: " << m_rclus);
96  ATH_MSG_INFO(" pT faction min: " << m_ptfrac);
97  ATH_MSG_INFO(" Jet builder: " << m_bld.name());
98 }
99 
100 //**********************************************************************
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
xAOD::JetAlgorithmType::kt_algorithm
@ kt_algorithm
Definition: JetContainerInfo.h:31
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
PseudoJetContainer
Definition: PseudoJetContainer.h:48
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
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
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
JetTrimmer::m_ptfrac
float m_ptfrac
Definition: JetTrimmer.h:48
JetTrimmer::groom
int groom(const xAOD::Jet &jin, const PseudoJetContainer &, xAOD::JetContainer &jout) const
Transform jet.
Definition: JetTrimmer.cxx:53
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
JetTrimmer::JetTrimmer
JetTrimmer(const std::string &name)
Definition: JetTrimmer.cxx:21
xAOD::Jet_v1::setAttribute
void setAttribute(const std::string &name, const T &v)
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:192
IJetPseudojetRetriever::pseudojet
virtual const fastjet::PseudoJet * pseudojet(const xAOD::Jet &jet) const =0
Retrieve the pseudojet associate with a jet.
xAOD::JetTransform::Trim
@ Trim
Definition: JetContainerInfo.h:115
xAOD::Jet_v1
Class describing a jet.
Definition: Jet_v1.h:57
JetTrimmer::print
void print() const
Print the state of the tool.
Definition: JetTrimmer.cxx:94
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
JetTrimmer.h
JetTrimmer::initialize
StatusCode initialize()
Dummy implementation of the initialisation function.
Definition: JetTrimmer.cxx:34
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
python.StandardJetMods.Filter
Filter
Definition: StandardJetMods.py:36
JetTrimmer::m_bld
ToolHandle< IJetFromPseudojet > m_bld
Definition: JetTrimmer.h:49
JetTrimmer::m_rclus
float m_rclus
Definition: JetTrimmer.h:47
JetTrimmer::~JetTrimmer
~JetTrimmer()