ATLAS Offline Software
JetSoftDrop.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 // JetSoftDrop.cxx
6 
7 #include "JetRec/JetSoftDrop.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 
16 using std::setw;
17 using fastjet::PseudoJet;
18 using xAOD::JetContainer;
19 
20 //**********************************************************************
21 
22 JetSoftDrop::JetSoftDrop(const std::string& name)
23  : AsgTool(name), m_bld("",this) {
24  //ZCut = 0.1, Beta = 0 are "CMS-like" parameters for SoftDrop:
25  //https://indico.cern.ch/event/439039/contributions/2223279/attachments/1311773/1963161/BOOST16_toptagging_CMS.pdf
26  declareProperty("ZCut", m_zcut =0.1);
27  declareProperty("Beta", m_beta =0.0);
28  declareProperty("R0", m_R0 =1.0);
29  declareProperty("JetBuilder", m_bld);
30 }
31 
32 //**********************************************************************
33 
34 JetSoftDrop::~JetSoftDrop() = default;
35 
36 //**********************************************************************
37 
39  if ( m_zcut < 0.0 || m_zcut > 10.0 ) {
40  ATH_MSG_ERROR("Invalid value for ZCut " << m_zcut);
41  return StatusCode::FAILURE;
42  }
43  if ( m_beta < 0.0 || m_beta > 10.0 ) {
44  ATH_MSG_ERROR("Invalid value for Beta " << m_beta);
45  return StatusCode::FAILURE;
46  }
47  if ( m_R0 < 0.0 || m_R0 > 10.0 ) {
48  ATH_MSG_ERROR("Invalid value for R0 " << m_R0);
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 
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 
74  //configure soft drop tool
75  //http://fastjet.hepforge.org/svn/contrib/contribs/RecursiveTools/tags/1.0.0/SoftDrop.hh
77  fastjet::contrib::SoftDrop softdropper(m_beta, m_zcut, m_R0);
78  PseudoJet pjsoftdrop = softdropper(*ppjin);
79  int npsoftdrop = pjsoftdrop.pieces().size();
80  xAOD::Jet* pjet = m_bld->add(pjsoftdrop, pjContainer, jets, &jin);
81  if ( pjet == nullptr ) {
82  ATH_MSG_ERROR("Unable to add jet to container");
83  return 1;
84  }
85  //pjet->SetAttribute<int>("TransformType", xAOD::JetTransform::SoftDrop); //xAOD::JetTransform::SoftDrop probably doesn't exist yet
86  pjet->setAttribute<float>("ZCut", m_zcut);
87  pjet->setAttribute<float>("SoftDropBeta", m_beta);
88  pjet->setAttribute<float>("SoftDropR0", m_R0);
89  pjet->setAttribute<int>("NSoftDropSubjets", npsoftdrop);
90 
91  ATH_MSG_DEBUG("Properties after softdrop:");
92  ATH_MSG_DEBUG(" ncon: " << pjsoftdrop.constituents().size() << "/"
93  << ppjin->constituents().size());
94  ATH_MSG_DEBUG(" nsub: " << npsoftdrop);
95  ATH_MSG_DEBUG("Added jet to container.");
96 
97  return 0;
98 }
99 
100 //**********************************************************************
101 
102 void JetSoftDrop::print() const {
103  ATH_MSG_INFO(" Asymmetry measure min: " << m_zcut);
104  ATH_MSG_INFO(" Angular exponent: " << m_beta);
105  ATH_MSG_INFO(" Characteristic jet radius: " << m_R0);
106  ATH_MSG_INFO(" Jet builder: " << m_bld.name());
107 }
108 
109 //**********************************************************************
JetSoftDrop::groom
int groom(const xAOD::Jet &jin, const PseudoJetContainer &, xAOD::JetContainer &jout) const
Transform jet.
Definition: JetSoftDrop.cxx:60
JetSoftDrop::m_bld
ToolHandle< IJetFromPseudojet > m_bld
Definition: JetSoftDrop.h:60
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
JetSoftDrop::~JetSoftDrop
~JetSoftDrop()
JetSoftDrop::JetSoftDrop
JetSoftDrop(const std::string &name)
Definition: JetSoftDrop.cxx:22
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
JetSoftDrop::m_zcut
float m_zcut
Definition: JetSoftDrop.h:57
PseudoJetContainer
Definition: PseudoJetContainer.h:48
JetSoftDrop::print
void print() const
Print the state of the tool.
Definition: JetSoftDrop.cxx:102
JetSoftDrop.h
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
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
JetSoftDrop::m_R0
float m_R0
Definition: JetSoftDrop.h:59
xAOD::Jet_v1::setAttribute
void setAttribute(const std::string &name, const T &v)
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
JetSoftDrop::initialize
StatusCode initialize()
Dummy implementation of the initialisation function.
Definition: JetSoftDrop.cxx:38
IJetPseudojetRetriever::pseudojet
virtual const fastjet::PseudoJet * pseudojet(const xAOD::Jet &jet) const =0
Retrieve the pseudojet associate with a jet.
xAOD::Jet_v1
Class describing a jet.
Definition: Jet_v1.h:57
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
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
JetSoftDrop::m_beta
float m_beta
Definition: JetSoftDrop.h:58