ATLAS Offline Software
Loading...
Searching...
No Matches
RunKLFitterAlg.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
7
8#ifndef KLFITTERNANALYSISALGORITHMS_RUNKLFITTERALG_H_
9#define KLFITTERNANALYSISALGORITHMS_RUNKLFITTERALG_H_
10
11// Algorithm includes
18
19// Framework includes
25
28
29// Externals
33#include "KLFitter/BoostedLikelihoodTopLeptonJets.h"
34#include "KLFitter/DetectorAtlas_8TeV.h"
35#include "KLFitter/Fitter.h"
36#include "KLFitter/LikelihoodTTHLeptonJets.h"
37#include "KLFitter/LikelihoodTTZTrilepton.h"
38#include "KLFitter/LikelihoodTopAllHadronic.h"
39#include "KLFitter/LikelihoodTopLeptonJets.h"
40#include "KLFitter/LikelihoodTopLeptonJets_Angular.h"
41#include "KLFitter/LikelihoodTopLeptonJets_JetAngles.h"
42#include "KLFitter/Permutations.h"
43
44namespace EventReco {
45
46class RunKLFitterAlg final : public EL::AnaAlgorithm {
47
48 public:
50 virtual StatusCode initialize() final;
51 virtual StatusCode execute() final;
52
53 private:
54 StatusCode execute_syst(const CP::SystematicSet &sys);
55 StatusCode add_leptons(
56 const std::vector<const xAOD::Electron *> &selected_electrons,
57 const std::vector<const xAOD::Muon *> &selected_muons,
58 KLFitter::Particles *myParticles);
59
60 StatusCode add_jets(const std::vector<const xAOD::Jet *> &selected_jets,
61 KLFitter::Particles *myParticles);
62
63 StatusCode setJetskLeadingN(const std::vector<const xAOD::Jet *> &jets,
64 KLFitter::Particles *inputParticles,
65 const size_t njets);
66
67 StatusCode retrieveEfficiencies(const xAOD::Jet *jet, float *eff,
68 float *ineff);
69
70 StatusCode setJetskBtagPriority(const std::vector<const xAOD::Jet *> &jets,
71 KLFitter::Particles *inputParticles,
72 const size_t maxJets);
73
74 StatusCode evaluatePermutations(const CP::SystematicSet &sys,
75 const std::vector<size_t> &electron_indices,
76 const std::vector<size_t> &muon_indices,
77 const std::vector<size_t> &jet_indices);
78
79 template <typename T>
80 std::vector<const T *> sortPt(const std::vector<const T *> &particles,
81 std::vector<size_t> &indices) {
82 std::vector<std::pair<const T *, size_t>> particle_index(particles.size());
83 size_t indx{0};
84 for (const T *const p : particles) {
85 particle_index[indx] = {p, indx};
86 ++indx;
87 }
89 particle_index.begin(), particle_index.end(),
90 [](std::pair<const T *, size_t> &x, std::pair<const T *, size_t> &y) {
91 return x.first->pt() > y.first->pt();
92 });
93 std::vector<const T *> sorted_particles(particles.size());
94 indx = 0;
95 indices.resize(particles.size());
96 for (auto &elem : particle_index) {
97 sorted_particles[indx] = elem.first;
98 indices[indx] = elem.second;
99 ++indx;
100 }
101 return sorted_particles;
102 }
103
104 // systematics
106
107 // inputs needed for reconstruction
109 this, "electrons", "", "the electron container to use"};
111 this, "electronSelection", "", "the selection on the input electrons"};
112
114 this, "muons", "", "the muon container to use"};
116 this, "muonSelection", "", "the selection on the input muons"};
117
119 this, "jets", "", "the jet container to use"};
121 "the selection on the input jets"};
122
124 this, "met", "", "the MET container to use"};
125
127 this, "eventInfo", "EventInfo",
128 "the EventInfo container to read selection deciosions from"};
129
130 // output container
133 m_outHandle{this, "result", "KLFitterResult_%SYS%",
134 "the output KLFitterResultContainer"};
135
136 CP::SysReadSelectionHandle m_selection{this, "selectionDecorationName", "",
137 "Name of the selection on which this "
138 "KLFitter instance is allowed to run"};
139
140 // configurable properties
141 Gaudi::Property<std::string> m_leptonType{this, "LeptonType", "kUndefined",
142 "Define the lepton type"};
143 Gaudi::Property<std::string> m_LHType{this, "LHType", "kUndefined",
144 "Define the Likelihood type"};
145 Gaudi::Property<std::string> m_transferFunctionsPath{
146 this, "TransferFunctionsPath",
147 "dev/AnalysisTop/KLFitterTFs/mc12a/akt4_LCtopo_PP6/",
148 "Path to transfer functions"};
149 Gaudi::Property<std::string> m_jetSelectionMode{
150 this, "JetSelectionMode", "kBtagPriorityFourJets",
151 "Define the behavior for selecting jets"};
152 Gaudi::Property<std::string> m_bTaggingMethod{
153 this, "BTaggingMethod", "kNotag",
154 "Method for accounting b-tagging information"};
155 Gaudi::Property<std::string> m_bTagDecoration{
156 this, "BTaggingDecoration", "",
157 "Name of the btag decision decoration for jets"};
158 Gaudi::Property<std::string> m_METterm{this, "METterm", "Final",
159 "Which MET term should be used"};
160 Gaudi::Property<float> m_massTop{
161 this, "TopMass", 172.5,
162 "The mass of top quark used in KLFitter likelihood (assuming the fixed "
163 "m_top mode is used)"};
164 Gaudi::Property<bool> m_fixedTopMass{
165 this, "TopMassFixed", true,
166 "If the top quark mass is fixed in the likelihood to the value of "
167 "TopMass parameter"};
168 Gaudi::Property<bool> m_saveAllPermutations{
169 this, "SaveAllPermutations", false,
170 "Whether to store only the permutation with highest KLFitter event "
171 "probability, or all"};
172 Gaudi::Property<bool> m_failOnLessThanXJets{
173 this, "FailOnLessThanXJets", false,
174 "Fail if kLeadingX or kBtagPriorityXJets is set and the number of jets "
175 "in the event is less than X"};
176
180
181 bool m_useBtagPriority{false};
183
184 std::unique_ptr<KLFitter::Fitter> m_myFitter;
185 std::unique_ptr<KLFitter::DetectorAtlas_8TeV> m_myDetector;
187 KLFitter::LikelihoodBase::BtaggingMethod m_bTaggingMethodEnum{};
188
189 KLFitter::LikelihoodTopLeptonJets::LeptonType m_leptonTypeKLFitterEnum{};
190 KLFitter::LikelihoodTTHLeptonJets::LeptonType m_leptonTypeKLFitterEnum_TTH{};
191 KLFitter::LikelihoodTopLeptonJets_JetAngles::LeptonType
193 KLFitter::LikelihoodTopLeptonJets_Angular::LeptonType
195 KLFitter::LikelihoodTTZTrilepton::LeptonType m_leptonTypeKLFitterEnum_TTZ{};
196 KLFitter::BoostedLikelihoodTopLeptonJets::LeptonType
198
199 std::unique_ptr<KLFitter::LikelihoodTopLeptonJets> m_myLikelihood;
200 std::unique_ptr<KLFitter::LikelihoodTTHLeptonJets> m_myLikelihood_TTH;
201 std::unique_ptr<KLFitter::LikelihoodTopLeptonJets_JetAngles>
203 std::unique_ptr<KLFitter::LikelihoodTopLeptonJets_Angular>
205 std::unique_ptr<KLFitter::LikelihoodTTZTrilepton> m_myLikelihood_TTZ;
206 std::unique_ptr<KLFitter::LikelihoodTopAllHadronic>
208 std::unique_ptr<KLFitter::BoostedLikelihoodTopLeptonJets>
210
211 ToolHandle<IBTaggingEfficiencyTool> m_btagging_eff_tool{
212 this, "btagEffTool", "", "the b-tagging efficiency tool"};
213
214 std::unique_ptr<SG::ConstAccessor<char>> m_bTagDecoAcc;
215};
216} // namespace EventReco
217
218#endif
#define y
#define x
a class managing the property to configure the list of systematics to process
a data handle for reading systematics varied input data
a data handle for reading systematically varied selection properties from objects
a data handle for writing systematics varied input data
Class to wrap a set of SystematicVariations.
the (new) base class for EventLoop algorithms
AnaAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
constructor with parameters
StatusCode evaluatePermutations(const CP::SystematicSet &sys, const std::vector< size_t > &electron_indices, const std::vector< size_t > &muon_indices, const std::vector< size_t > &jet_indices)
CP::SysWriteHandle< xAOD::KLFitterResultContainer, xAOD::KLFitterResultAuxContainer > m_outHandle
KLFitter::LikelihoodTTZTrilepton::LeptonType m_leptonTypeKLFitterEnum_TTZ
virtual StatusCode initialize() final
StatusCode retrieveEfficiencies(const xAOD::Jet *jet, float *eff, float *ineff)
std::unique_ptr< KLFitter::LikelihoodTopLeptonJets_JetAngles > m_myLikelihood_JetAngles
CP::SysReadHandle< xAOD::EventInfo > m_eventInfoHandle
std::unique_ptr< KLFitter::LikelihoodTTHLeptonJets > m_myLikelihood_TTH
Gaudi::Property< std::string > m_transferFunctionsPath
StatusCode setJetskLeadingN(const std::vector< const xAOD::Jet * > &jets, KLFitter::Particles *inputParticles, const size_t njets)
CP::SysListHandle m_systematicsList
std::unique_ptr< SG::ConstAccessor< char > > m_bTagDecoAcc
std::unique_ptr< KLFitter::DetectorAtlas_8TeV > m_myDetector
CP::SysReadSelectionHandle m_muonSelection
StatusCode execute_syst(const CP::SystematicSet &sys)
CP::SysReadHandle< xAOD::ElectronContainer > m_electronsHandle
KLFEnums::JetSelectionMode m_jetSelectionModeEnum
std::unique_ptr< KLFitter::LikelihoodTTZTrilepton > m_myLikelihood_TTZ
KLFEnums::LeptonType m_leptonTypeEnum
std::vector< const T * > sortPt(const std::vector< const T * > &particles, std::vector< size_t > &indices)
Gaudi::Property< bool > m_fixedTopMass
CP::SysReadHandle< xAOD::JetContainer > m_jetsHandle
std::unique_ptr< KLFitter::LikelihoodTopLeptonJets > m_myLikelihood
Gaudi::Property< std::string > m_leptonType
KLFitter::BoostedLikelihoodTopLeptonJets::LeptonType m_leptonTypeKLFitterEnum_BoostedLJets
KLFitter::LikelihoodTopLeptonJets_Angular::LeptonType m_leptonTypeKLFitterEnum_Angular
std::unique_ptr< KLFitter::BoostedLikelihoodTopLeptonJets > m_myLikelihood_BoostedLJets
Gaudi::Property< std::string > m_bTaggingMethod
StatusCode add_leptons(const std::vector< const xAOD::Electron * > &selected_electrons, const std::vector< const xAOD::Muon * > &selected_muons, KLFitter::Particles *myParticles)
Gaudi::Property< bool > m_saveAllPermutations
Gaudi::Property< std::string > m_METterm
CP::SysReadHandle< xAOD::MissingETContainer > m_metHandle
KLFEnums::Likelihood m_LHTypeEnum
CP::SysReadSelectionHandle m_selection
KLFitter::LikelihoodBase::BtaggingMethod m_bTaggingMethodEnum
ToolHandle< IBTaggingEfficiencyTool > m_btagging_eff_tool
Gaudi::Property< std::string > m_bTagDecoration
KLFitter::LikelihoodTopLeptonJets_JetAngles::LeptonType m_leptonTypeKLFitterEnum_JetAngles
StatusCode setJetskBtagPriority(const std::vector< const xAOD::Jet * > &jets, KLFitter::Particles *inputParticles, const size_t maxJets)
KLFEnums::JetSelectionMode m_jetSelectionModeKLFitterEnum
CP::SysReadHandle< xAOD::MuonContainer > m_muonsHandle
std::unique_ptr< KLFitter::LikelihoodTopLeptonJets_Angular > m_myLikelihood_Angular
KLFitter::LikelihoodTTHLeptonJets::LeptonType m_leptonTypeKLFitterEnum_TTH
Gaudi::Property< bool > m_failOnLessThanXJets
CP::SysReadSelectionHandle m_electronSelection
virtual StatusCode execute() final
Gaudi::Property< std::string > m_jetSelectionMode
KLFitter::LikelihoodTopLeptonJets::LeptonType m_leptonTypeKLFitterEnum
StatusCode add_jets(const std::vector< const xAOD::Jet * > &selected_jets, KLFitter::Particles *myParticles)
CP::SysReadSelectionHandle m_jetSelection
Gaudi::Property< std::string > m_LHType
std::unique_ptr< KLFitter::Fitter > m_myFitter
Gaudi::Property< float > m_massTop
std::unique_ptr< KLFitter::LikelihoodTopAllHadronic > m_myLikelihood_AllHadronic
Auxiliary container for xAOD::KLFitterResultContainer.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
Jet_v1 Jet
Definition of the current "jet version".
DataVector< xAOD::KLFitterResult > KLFitterResultContainer
Definition of the KLFitterResultContainer type.