ATLAS Offline Software
Loading...
Searching...
No Matches
MHTPufitFex.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3 */
4
5/******************************************************************************
6 * @package Trigger/TrigAlgorithms/TrigEFMissingET
7 * @class MHTPufitFex
8 *
9 * @brief Fex class for the mhtpufit algorithm, using the pufit procedure to
10 * correct the MET produced by a sum over jets
11 * @author Jon Burr
12 *****************************************************************************/
13
14#ifndef TRIGEFMISSINGET_MHTPUFITFEX_H
15#define TRIGEFMISSINGET_MHTPUFITFEX_H
16
17#include "FexBase.h"
23#include "Gaudi/Property.h"
24#include "GaudiKernel/SystemOfUnits.h"
25
26namespace HLT
27{
28 namespace MET
29 {
30 /**************************************************************************
31 * @class MHTPufitFex
32 *
33 * Class to calculate MET from a sum over jets using a correction derived
34 * with the pufit technique
35 *************************************************************************/
36 class MHTPufitFex : public FexBase
37 {
38 public:
40 MHTPufitFex(const std::string &name, ISvcLocator *pSvcLocator);
41
43 virtual StatusCode initialize() override;
44
45 private:
46 /************************************************************************
47 * Properties
48 ***********************************************************************/
51 this, "InputJetsName", "", "The input jet container"};
52
54 this, "InputJvtName", "Jvt", "The input JVT name"};
55
57 this, "InputName", "", "The input clusters or PFOs"};
58 Gaudi::Property<bool> m_jetCalibIncludesAreaSub{
59 this, "JetCalibIncludesAreaSub", true,
60 "Whether the calibration applied to the jets includes area subtraction"};
62 this, "JetEventShapeName", "",
63 "The name of the event shape container for the area correction"};
64
65 Gaudi::Property<float> m_nSigma{
66 this, "NSigma", 5, "Set the threshold at mean + NSigma*variance"};
67
68 Gaudi::Property<float> m_maxEta{
69 this, "MaxEta", 5, "The maximum eta range"};
70
71 Gaudi::Property<std::size_t> m_nEtaBins{
72 this, "NEtaBins", 14, "The number of eta bins"};
73 Gaudi::Property<std::size_t> m_nPhiBins{
74 this, "NPhiBins", 8, "The number of phi bins"};
75
76 Gaudi::Property<float> m_trimFraction{
77 this, "TrimFraction", 0.9,
78 "The fraction of bins to use when calculating the mean and variance"};
79
80 Gaudi::Property<float> m_caloNoise{
81 this, "CaloNoise", 50,
82 "The coefficient of the noise term in the calorimeter resolution estimate [MeV]"};
83
84 Gaudi::Property<float> m_caloStoch{
85 this, "CaloStochastic", 15.81,
86 "The coefficient of the stochastic term in the calorimeter resolution estimate [MeV^1/2]"};
87
88 Gaudi::Property<float> m_constraintImportance{
89 this, "ConstraintImportance", 1,
90 "The relative importance of the two constraints to the fit"};
91 // Properties relating to jet selection
92 Gaudi::Property<float> m_jvtCut{
93 this, "JvtCut", 0.59, "The JVT selection in the central region"};
94 Gaudi::Property<float> m_minPt{
95 this, "MinPt", 20*Gaudi::Units::GeV, "The minimum pT (in the central region)"};
96 Gaudi::Property<float> m_maxPt{
97 this, "MaxPt", 120*Gaudi::Units::GeV,
98 "The maximum pT (in the central region), above which the JVT selection is not applied"};
99 Gaudi::Property<float> m_forwardPt{
100 this, "ForwardPt", 30*Gaudi::Units::GeV, "The minimum pt in the forward region"};
101 Gaudi::Property<float> m_centralEta{
102 this, "CentralEta", 2.4, "The boundary between the central and border regions"};
103 Gaudi::Property<bool> m_useDetectorEta{
104 this, "UseDetectorEta", true, "Whether to use the 'DectectorEta' value to select central/forward jets"};
105 /************************************************************************
106 * Internal functions
107 ***********************************************************************/
114 virtual StatusCode fillMET(
116 const EventContext &context,
117 MonGroupBuilder &monitors) const override;
118
123 double getSigma(const SignedKinematics &kin) const;
124 }; //> end class MHTPufitFex
125 } // namespace MET
126} // namespace HLT
127
128#endif //> !TRIGEFMISSINGET_MHTPUFITFEX_H
Property holding a SG store/key/clid/attr name from which a ReadDecorHandle is made.
FexBase(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
Definition FexBase.cxx:43
Gaudi::Property< float > m_minPt
Definition MHTPufitFex.h:94
virtual StatusCode initialize() override
Initialize the fex.
Gaudi::Property< bool > m_useDetectorEta
virtual StatusCode fillMET(xAOD::TrigMissingET &met, const EventContext &context, MonGroupBuilder &monitors) const override
Calculate and fill the output MET value.
Gaudi::Property< bool > m_jetCalibIncludesAreaSub
Definition MHTPufitFex.h:58
Gaudi::Property< float > m_trimFraction
The trimming fraction.
Definition MHTPufitFex.h:76
Gaudi::Property< float > m_nSigma
The sigma threshold.
Definition MHTPufitFex.h:65
Gaudi::Property< float > m_caloStoch
The coefficient of the stochastic term in the calo resolution estimate.
Definition MHTPufitFex.h:84
Gaudi::Property< float > m_forwardPt
Definition MHTPufitFex.h:99
Gaudi::Property< float > m_centralEta
Gaudi::Property< float > m_caloNoise
The coefficient of the noise term in the calo resolution estimate.
Definition MHTPufitFex.h:80
Gaudi::Property< float > m_constraintImportance
The relative constraint importance.
Definition MHTPufitFex.h:88
SG::ReadDecorHandleKey< xAOD::JetContainer > m_inputJvtKey
The input JVT decoration.
Definition MHTPufitFex.h:53
SG::ReadHandleKey< xAOD::EventShape > m_rhoKey
Definition MHTPufitFex.h:61
SG::ReadHandleKey< xAOD::IParticleContainer > m_inputKey
The input clusters or PFOs.
Definition MHTPufitFex.h:56
Gaudi::Property< float > m_maxEta
The eta range of the grid.
Definition MHTPufitFex.h:68
Gaudi::Property< std::size_t > m_nPhiBins
Definition MHTPufitFex.h:73
Gaudi::Property< std::size_t > m_nEtaBins
The number of bins in eta.
Definition MHTPufitFex.h:71
Gaudi::Property< float > m_jvtCut
Definition MHTPufitFex.h:92
double getSigma(const SignedKinematics &kin) const
Calculate the estimate on the variance of a tower.
MHTPufitFex(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
SG::ReadHandleKey< xAOD::JetContainer > m_inputJetsKey
Input objects.
Definition MHTPufitFex.h:50
Gaudi::Property< float > m_maxPt
Definition MHTPufitFex.h:96
Class to describe the kinematics of an object that can have negative energies.
Property holding a SG store/key/clid/attr name from which a ReadDecorHandle is made.
Property holding a SG store/key/clid from which a ReadHandle is made.
It used to be useful piece of code for replacing actual SG with other store of similar functionality ...
TrigMissingET_v1 TrigMissingET
Define the most recent version of the TrigMissingET class.