ATLAS Offline Software
Loading...
Searching...
No Matches
JetRecursiveSoftDrop.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023c CERN for the benefit of the ATLAS collaboration
3*/
4
5// JetRecursiveSoftDrop.cxx
7#include <iomanip>
8#include "fastjet/PseudoJet.hh"
9#include "fastjet/JetDefinition.hh"
11
12
13using std::setw;
14using fastjet::PseudoJet;
16
17//**********************************************************************
18
20 : AsgTool(name), m_bld("",this) {
21 declareProperty("ZCut", m_zcut =0.1);
22 declareProperty("Beta", m_beta =0.0);
23 declareProperty("N", m_N =-1); //default to infinite layers
24 declareProperty("R0", m_R0 =1.0);
25 declareProperty("JetBuilder", m_bld);
26}
27
28//**********************************************************************
29
31
32//**********************************************************************
33
35 if ( m_zcut < 0.0 || m_zcut > 10.0 ) {
36 ATH_MSG_ERROR("Invalid value for ZCut " << m_zcut);
37 return StatusCode::FAILURE;
38 }
39 if ( m_beta < 0.0 || m_beta > 10.0 ) {
40 ATH_MSG_ERROR("Invalid value for Beta " << m_beta);
41 return StatusCode::FAILURE;
42 }
43 if ( m_N < -1.0) {
44 ATH_MSG_ERROR("Invalid value for N " << m_N);
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 recursive soft drop tool
75 //https://fastjet.hepforge.org/trac/browser/contrib/contribs/RecursiveTools/tags/2.0.0-beta1
77 fastjet::contrib::RecursiveSoftDrop softdropper(m_beta, m_zcut, m_N, 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<float>("RSoftDropZCut", m_zcut);
86 pjet->setAttribute<float>("RSoftDropBeta", m_beta);
87 pjet->setAttribute<int>("RSoftDropN", m_N);
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
96 ATH_MSG_DEBUG("Added jet to container.");
97
98 return 0;
99}
100
101//**********************************************************************
102
104 ATH_MSG_INFO(" Asymmetry measure min: " << m_zcut);
105 ATH_MSG_INFO(" Angular exponent: " << m_beta);
106 ATH_MSG_INFO(" Recursive layers: " << m_N);
107 ATH_MSG_INFO(" Characteristic jet radius: " << m_R0);
108 ATH_MSG_INFO(" Jet builder: " << m_bld.name());
109}
110
111//**********************************************************************
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
virtual const IJetPseudojetRetriever * pseudojetRetriever() const
Return the pseudojet retriever associated with this tool.
virtual const fastjet::PseudoJet * pseudojet(const xAOD::Jet &jet) const =0
Retrieve the pseudojet associate with a jet.
ToolHandle< IJetFromPseudojet > m_bld
JetRecursiveSoftDrop(const std::string &name)
StatusCode initialize()
Dummy implementation of the initialisation function.
void print() const
Print the state of the tool.
int groom(const xAOD::Jet &jin, const PseudoJetContainer &, xAOD::JetContainer &jout) const
Transform jet.
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
void setAttribute(const std::string &name, const T &v)
Jet_v1 Jet
Definition of the current "jet version".
JetContainer_v1 JetContainer
Definition of the current "jet container version".