ATLAS Offline Software
Loading...
Searching...
No Matches
InSituCalibStep.h
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
5*/
6
7// InSituCalibStep.h
8// Header file for class InSituCalibStep
9// Author: Fabrice Balli <fabrice.balli@cern.ch>
11#ifndef JETCALIBTOOLS_INSITUCALIBSTEP_H
12#define JETCALIBTOOLS_INSITUCALIBSTEP_H 1
13
14#include <string.h>
15
16#include <TString.h>
17#include <TEnv.h>
18
19#include "AsgTools/AsgTool.h"
21#include "AsgTools/ToolHandle.h"
24
26
32
34 : public asg::AsgTool,
35 virtual public IJetCalibStep {
37
38public:
40 InSituCalibStep(const std::string& name = "InSituCalibStep");
41
42 virtual StatusCode initialize() override;
43 virtual StatusCode calibrate(xAOD::JetContainer&) const override;
44
45
46private:
47
48 Gaudi::Property<bool> m_CalibrateMC {this, "CalibrateMC", false, "force Insitu step for MC sample"};
49 Gaudi::Property<bool> m_isMC {this, "isMC", false, "isMC"};
50 Gaudi::Property<bool> m_useOriginalHistCombination {this, "useOriginalHistCombination",false,"Use original method to combine histograms."};
51
52 Gaudi::Property<std::string> m_jetInScale {this, "InScale", "JetGSCScaleMomentum", "Starting jet scale"};
53 Gaudi::Property<std::string> m_jetOutScale {this, "OutScale", "JetInsituScaleMomentum", "Ending jet scale"};
54 // retrieves in situ correction
55 StatusCode getInsituCorr(const xAOD::Jet& jet, JetHelper::JetContext& jc, unsigned int periodIndex, double &scale) const;
56 // Relative calibration (derived with eta intercalibration)
57 //ToolHandle<JetHelper::IVarTool> m_histTool_EtaInter{this, "HistoReaderEtaInter", "HistoInput2D", "Instance of HistoInput2D for reading histogram"};
58 ToolHandleArray<JetHelper::IVarTool> m_histTool_EtaInter{this, "HistoReaderEtaInter", {}, "Instance of HistoInput2D for reading histogram"};
59 // Absolute calibration (derived with |eta| < 0.8)
60 //ToolHandle<JetHelper::IVarTool> m_histTool_Abs{this, "HistoReaderAbs", "HistoInput1D", "Instance of HistoInput1D for reading histogram"};
61 ToolHandleArray<JetHelper::IVarTool> m_histTool_Abs{this, "HistoReaderAbs", {}, "Instance of HistoInput1D for reading histogram"};
62 // needed for the combined histogram: pT
63 ToolHandle<JetHelper::IVarTool> m_vartool1 {this, "vartool1", "VarTool", "InputVariable instance" };
64 // needed for the combined histogram: eta
65 ToolHandle<JetHelper::IVarTool> m_vartool2 {this, "vartool2", "VarTool", "InputVariable instance" };
66 // vector of run numbers
67 Gaudi::Property<std::vector<unsigned int> > m_RunNumBoundaries{this, "RunNumbers", {} ,""};
68 // combine Relative and Absolute calibrations
69 std::unique_ptr<const TH2> combineCalibration(const TH2* h2d, const TH1* h);
70 // combined in situ correction histogram
71 std::vector< std::unique_ptr<const TH2> > m_insituCorr_vec;
72 // maximum eta of combined histogram
73 std::vector<double> m_etaMax_vec;
74 // minimum eta of combined histogram
75 std::vector<double> m_etaMin_vec;
76 // maximum pt of combined histogram
77 std::vector<double> m_ptMax_vec;
78 // minimum pt of combined histogram
79 std::vector<double> m_ptMin_vec;
80 //ReadHandleKey for event info
81 SG::ReadHandleKey<xAOD::EventInfo> m_evtInfoKey{this, "EventInfoKey", "EventInfo"};
82 // Method to get runNumber
83 StatusCode retrieveEventInfo(unsigned int &r) const;
84
85
86};
87#endif
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
Property holding a SG store/key/clid/attr name from which a ReadDecorHandle is made.
Header file for AthHistogramAlgorithm.
std::unique_ptr< const TH2 > combineCalibration(const TH2 *h2d, const TH1 *h)
ToolHandleArray< JetHelper::IVarTool > m_histTool_Abs
Gaudi::Property< std::string > m_jetInScale
Gaudi::Property< bool > m_CalibrateMC
StatusCode retrieveEventInfo(unsigned int &r) const
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
virtual StatusCode calibrate(xAOD::JetContainer &) const override
Apply calibration to a jet container.
std::vector< double > m_ptMax_vec
std::vector< double > m_ptMin_vec
std::vector< std::unique_ptr< const TH2 > > m_insituCorr_vec
Gaudi::Property< std::string > m_jetOutScale
std::vector< double > m_etaMin_vec
ToolHandleArray< JetHelper::IVarTool > m_histTool_EtaInter
ToolHandle< JetHelper::IVarTool > m_vartool2
StatusCode getInsituCorr(const xAOD::Jet &jet, JetHelper::JetContext &jc, unsigned int periodIndex, double &scale) const
std::vector< double > m_etaMax_vec
SG::ReadHandleKey< xAOD::EventInfo > m_evtInfoKey
ToolHandle< JetHelper::IVarTool > m_vartool1
InSituCalibStep(const std::string &name="InSituCalibStep")
Constructor with parameters:
Gaudi::Property< bool > m_useOriginalHistCombination
Gaudi::Property< std::vector< unsigned int > > m_RunNumBoundaries
Gaudi::Property< bool > m_isMC
Class JetContext Designed to read AOD information related to the event, N vertices,...
Definition JetContext.h:24
Property holding a SG store/key/clid from which a ReadHandle is made.
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
int r
Definition globals.cxx:22
Jet_v1 Jet
Definition of the current "jet version".
JetContainer_v1 JetContainer
Definition of the current "jet container version".