ATLAS Offline Software
Loading...
Searching...
No Matches
JMSCalibStep.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
7
8JMSCalibStep::JMSCalibStep(const std::string& name)
9 : asg::AsgTool( name ){ }
10
12
13 ATH_MSG_INFO("Initializing the JMS Calibration tool");
14
15 ATH_MSG_DEBUG("Reading from " << m_jetInScale << " and writing to " << m_jetOutScale);
16
17 ATH_CHECK(m_histTool.retrieve());
18
19 // Determine the maximum value on the z-azis for which the calibration will be applied, this should be eta
20 JetHelper::HistoInputBase* histoTool = dynamic_cast<JetHelper::HistoInputBase*>( &(*m_histTool) );
21 if(histoTool){
22 TH1 *h = &histoTool->getHistogram();
23 m_maxEta = h->GetZaxis()->GetBinLowEdge(h->GetNbinsZ()+1);
24 }
25 else{
26 ATH_MSG_ERROR ("Dynamic cast to JetHelper::HistoInputBase failed in JMSCalibStep!");
27 }
28
29 return StatusCode::SUCCESS;
30
31}
32
34
35 ATH_MSG_DEBUG("Calibrating jet mass");
36
38
39 for(xAOD::Jet* jet: jets){
40
41 const xAOD::JetFourMom_t jetStartP4 = jet->getAttribute<xAOD::JetFourMom_t>(m_jetInScale);
42 jet->setJetP4(jetStartP4);
43
44 // Get the calibration factor (only for jets if the pT or energy is above threshold and eta within the histogram z-axis range)
45 double massFactor = 1.0;
47 if(m_varToolX->getValue(*jet,jc) >= m_minValue_JMS){
48 if(m_varToolZ->getValue(*jet,jc) <= m_maxEta){
49 massFactor = m_histTool->getValue(*jet, jc);
50 }
51 }
52
53 // Calculate the corrected mass
54 double mass_corr = jetStartP4.mass();
55
56 if(massFactor != 0){
57 mass_corr = jetStartP4.mass()/massFactor;
58 }
59 // Protection for very large masses
60 if (mass_corr > jet->e()){
61 mass_corr = jet->m();
62 }
63
64 double pT_corr = jetStartP4.pt();
65 // For small-R jet mass calibrations, keep the pT value fixed
66 if(!m_pTfixed){
67 pT_corr = std::sqrt(jetStartP4.e()*jetStartP4.e()-mass_corr*mass_corr)/std::cosh( jetStartP4.eta() );
68 }
69
70 // Set the four-vector to the calibrated values
71 xAOD::JetFourMom_t calibP4 = xAOD::JetFourMom_t(pT_corr, jetStartP4.eta(), jetStartP4.phi(), mass_corr);
72 jmsScaleMomAcc.setAttribute(*jet, calibP4);
73 jet->setJetP4(calibP4);
74 }
75
76 return StatusCode::SUCCESS;
77
78}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
Header file for AthHistogramAlgorithm.
ToolHandle< JetHelper::IVarTool > m_varToolZ
Variable used to enforce eta threshold.
ToolHandle< JetHelper::IVarTool > m_histTool
3D histogram containing mass calibration values
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< std::string > m_jetInScale
Gaudi::Property< bool > m_pTfixed
Gaudi::Property< std::string > m_jetOutScale
JMSCalibStep(const std::string &name="JMSCalibStep")
ToolHandle< JetHelper::IVarTool > m_varToolX
Variable used to enforce minimum pT or energy for calibration.
virtual StatusCode calibrate(xAOD::JetContainer &) const override
Apply calibration to a jet container.
Gaudi::Property< float > m_minValue_JMS
Properties.
Class HistoInputBase This class implement common function used by HistoInput1D and HistoInput2D.
TH1 & getHistogram()
Returns the underlying histogram.
Class JetContext Designed to read AOD information related to the event, N vertices,...
Definition JetContext.h:24
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
void setAttribute(SG::AuxElement &p, const TYPE &v) const
Jet_v1 Jet
Definition of the current "jet version".
JetContainer_v1 JetContainer
Definition of the current "jet container version".
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17