CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Types | Private Attributes
METCorrectionAlgorithm Class Reference

#include <METCorrectionAlgorithm.h>

Classes

struct  type2BinningEntryType
 

Public Member Functions

CorrMETData compMETCorrection (edm::Event &)
 
 METCorrectionAlgorithm (const edm::ParameterSet &, edm::ConsumesCollector &&iConsumesCollector)
 
 ~METCorrectionAlgorithm ()
 

Private Types

typedef std::vector< edm::InputTagvInputTag
 
typedef std::vector< std::string > vstring
 

Private Attributes

bool applyType0Corrections_
 
bool applyType1Corrections_
 
bool applyType2Corrections_
 
std::vector< edm::EDGetTokenT< CorrMETData > > chsSumTokens_
 
double type0Cuncl_
 
double type0Rsoft_
 
std::vector< edm::EDGetTokenT< CorrMETData > > type1Tokens_
 
std::vector< type2BinningEntryType * > type2Binning_
 

Detailed Description

Algorithm for o propagating jet energy corrections to MET (Type 1 MET corrections) o calibrating momentum of particles not within jets ("unclustered energy") and propagating those corrections to MET (Type 2 MET corrections)

Authors
Michael Schmitt, Richard Cavanaugh, The University of Florida Florent Lacroix, University of Illinois at Chicago Christian Veelken, LLR Tai Sakuma, Texas A&M

Definition at line 32 of file METCorrectionAlgorithm.h.

Member Typedef Documentation

◆ vInputTag

typedef std::vector<edm::InputTag> METCorrectionAlgorithm::vInputTag
private

Definition at line 40 of file METCorrectionAlgorithm.h.

◆ vstring

typedef std::vector<std::string> METCorrectionAlgorithm::vstring
private

Definition at line 52 of file METCorrectionAlgorithm.h.

Constructor & Destructor Documentation

◆ METCorrectionAlgorithm()

METCorrectionAlgorithm::METCorrectionAlgorithm ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iConsumesCollector 
)
explicit

Definition at line 9 of file METCorrectionAlgorithm.cc.

References applyType0Corrections_, applyType1Corrections_, applyType2Corrections_, looper::cfg, chsSumTokens_, SimL1EmulatorRepack_Full_cff::inputTag, patCaloMETCorrections_cff::srcType1Corrections, correctionTermsCaloMet_cff::srcUnclEnergySums, AlCaHLTBitMon_QueryRunRegistry::string, type0Cuncl_, type0Rsoft_, type1Tokens_, type2Binning_, correctionTermsCaloMet_cff::type2CorrFormula, and correctionTermsCaloMet_cff::type2CorrParameter.

10  {
11  applyType1Corrections_ = cfg.getParameter<bool>("applyType1Corrections");
13  vInputTag srcType1Corrections = cfg.getParameter<vInputTag>("srcType1Corrections");
14  for (vInputTag::const_iterator inputTag = srcType1Corrections.begin(); inputTag != srcType1Corrections.end();
15  ++inputTag) {
16  type1Tokens_.push_back(iConsumesCollector.consumes<CorrMETData>(*inputTag));
17  }
18  }
19 
20  applyType2Corrections_ = cfg.getParameter<bool>("applyType2Corrections");
22  vInputTag srcUnclEnergySums = cfg.getParameter<vInputTag>("srcUnclEnergySums");
23 
24  if (cfg.exists("type2Binning")) {
25  typedef std::vector<edm::ParameterSet> vParameterSet;
26  vParameterSet cfgType2Binning = cfg.getParameter<vParameterSet>("type2Binning");
27  for (vParameterSet::const_iterator cfgType2BinningEntry = cfgType2Binning.begin();
28  cfgType2BinningEntry != cfgType2Binning.end();
29  ++cfgType2BinningEntry) {
30  type2Binning_.push_back(
31  new type2BinningEntryType(*cfgType2BinningEntry, srcUnclEnergySums, iConsumesCollector));
32  }
33  } else {
34  std::string type2CorrFormula = cfg.getParameter<std::string>("type2CorrFormula");
35  edm::ParameterSet type2CorrParameter = cfg.getParameter<edm::ParameterSet>("type2CorrParameter");
36  type2Binning_.push_back(
37  new type2BinningEntryType(type2CorrFormula, type2CorrParameter, srcUnclEnergySums, iConsumesCollector));
38  }
39  }
40 
42  cfg.exists("applyType0Corrections") ? cfg.getParameter<bool>("applyType0Corrections") : false;
44  vInputTag srcCHSSums = cfg.getParameter<vInputTag>("srcCHSSums");
45  for (vInputTag::const_iterator inputTag = srcCHSSums.begin(); inputTag != srcCHSSums.end(); ++inputTag) {
46  chsSumTokens_.push_back(iConsumesCollector.consumes<CorrMETData>(*inputTag));
47  }
48 
49  type0Rsoft_ = cfg.getParameter<double>("type0Rsoft");
50  type0Cuncl_ = 1.0;
52  if (cfg.exists("type2Binning"))
53  throw cms::Exception("Invalid Arg")
54  << "Currently, applyType0Corrections and type2Binning cannot be used together!";
55  std::string type2CorrFormula = cfg.getParameter<std::string>("type2CorrFormula");
56  if (!(type2CorrFormula == "A"))
57  throw cms::Exception("Invalid Arg") << "type2CorrFormula must be \"A\" if applyType0Corrections!";
58  edm::ParameterSet type2CorrParameter = cfg.getParameter<edm::ParameterSet>("type2CorrParameter");
59  type0Cuncl_ = type2CorrParameter.getParameter<double>("A");
60  }
61  }
62 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< edm::InputTag > vInputTag
std::vector< edm::EDGetTokenT< CorrMETData > > chsSumTokens_
std::vector< type2BinningEntryType * > type2Binning_
std::vector< edm::EDGetTokenT< CorrMETData > > type1Tokens_
a MET correction term
Definition: CorrMETData.h:14

◆ ~METCorrectionAlgorithm()

METCorrectionAlgorithm::~METCorrectionAlgorithm ( )

Definition at line 64 of file METCorrectionAlgorithm.cc.

References type2Binning_.

64  {
65  for (std::vector<type2BinningEntryType*>::const_iterator it = type2Binning_.begin(); it != type2Binning_.end();
66  ++it) {
67  delete (*it);
68  }
69 }
std::vector< type2BinningEntryType * > type2Binning_

Member Function Documentation

◆ compMETCorrection()

CorrMETData METCorrectionAlgorithm::compMETCorrection ( edm::Event evt)

Definition at line 71 of file METCorrectionAlgorithm.cc.

References applyType0Corrections_, applyType1Corrections_, applyType2Corrections_, chsSumTokens_, edm::Event::getByToken(), CorrMETData::mex, CorrMETData::mey, mathSSE::sqrt(), CorrMETData::sumet, type0Cuncl_, type0Rsoft_, type1Tokens_, type2Binning_, and trackerHitRTTI::vector.

71  {
72  CorrMETData metCorr;
73  metCorr.mex = 0.;
74  metCorr.mey = 0.;
75  metCorr.sumet = 0.;
76 
78  //--- sum all Type 0 MET correction terms
80  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = chsSumTokens_.begin();
81  corrToken != chsSumTokens_.end();
82  ++corrToken) {
83  evt.getByToken(*corrToken, chsSum);
84 
85  metCorr.mex += type0Cuncl_ * (1 - type0Rsoft_) * chsSum->mex;
86  metCorr.mey += type0Cuncl_ * (1 - type0Rsoft_) * chsSum->mey;
87  metCorr.sumet += type0Cuncl_ * (1 - type0Rsoft_) * chsSum->sumet;
88  }
89  }
90 
92  //--- sum all Type 1 MET correction terms
93  edm::Handle<CorrMETData> type1Correction;
94  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = type1Tokens_.begin();
95  corrToken != type1Tokens_.end();
96  ++corrToken) {
97  evt.getByToken(*corrToken, type1Correction);
98 
99  metCorr.mex += type1Correction->mex;
100  metCorr.mey += type1Correction->mey;
101  metCorr.sumet += type1Correction->sumet;
102  }
103  }
104 
106  //--- compute momentum sum of all "unclustered energy" in the event
107  //
108  // NOTE: calibration factors/formulas for Type 2 MET correction may depend on eta
109  // (like the jet energy correction factors do)
110  //
111 
112  for (std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
113  type2BinningEntry != type2Binning_.end();
114  ++type2BinningEntry) {
115  CorrMETData unclEnergySum;
116  edm::Handle<CorrMETData> unclEnergySummand;
117  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken =
118  (*type2BinningEntry)->corrTokens_.begin();
119  corrToken != (*type2BinningEntry)->corrTokens_.end();
120  ++corrToken) {
121  evt.getByToken(*corrToken, unclEnergySummand);
122 
123  unclEnergySum.mex += unclEnergySummand->mex;
124  unclEnergySum.mey += unclEnergySummand->mey;
125  unclEnergySum.sumet += unclEnergySummand->sumet;
126  }
127 
128  //--- calibrate "unclustered energy"
129  double unclEnergySumPt = sqrt(unclEnergySum.mex * unclEnergySum.mex + unclEnergySum.mey * unclEnergySum.mey);
130  double unclEnergyScaleFactor = (*type2BinningEntry)->binCorrFormula_->Eval(unclEnergySumPt);
131 
132  //--- MET balances momentum of reconstructed particles,
133  // hence correction to "unclustered energy" and corresponding Type 2 MET correction are of opposite sign
134  metCorr.mex -= (unclEnergyScaleFactor - 1.) * unclEnergySum.mex;
135  metCorr.mey -= (unclEnergyScaleFactor - 1.) * unclEnergySum.mey;
136  metCorr.sumet += (unclEnergyScaleFactor - 1.) * unclEnergySum.sumet;
137  }
138  }
139 
140  return metCorr;
141 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:536
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< edm::EDGetTokenT< CorrMETData > > chsSumTokens_
double sumet
Definition: CorrMETData.h:18
std::vector< type2BinningEntryType * > type2Binning_
std::vector< edm::EDGetTokenT< CorrMETData > > type1Tokens_
a MET correction term
Definition: CorrMETData.h:14
double mey
Definition: CorrMETData.h:16
double mex
Definition: CorrMETData.h:15

Member Data Documentation

◆ applyType0Corrections_

bool METCorrectionAlgorithm::applyType0Corrections_
private

Definition at line 45 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ applyType1Corrections_

bool METCorrectionAlgorithm::applyType1Corrections_
private

Definition at line 46 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ applyType2Corrections_

bool METCorrectionAlgorithm::applyType2Corrections_
private

Definition at line 47 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ chsSumTokens_

std::vector<edm::EDGetTokenT<CorrMETData> > METCorrectionAlgorithm::chsSumTokens_
private

Definition at line 42 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ type0Cuncl_

double METCorrectionAlgorithm::type0Cuncl_
private

Definition at line 50 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ type0Rsoft_

double METCorrectionAlgorithm::type0Rsoft_
private

Definition at line 49 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ type1Tokens_

std::vector<edm::EDGetTokenT<CorrMETData> > METCorrectionAlgorithm::type1Tokens_
private

Definition at line 43 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

◆ type2Binning_

std::vector<type2BinningEntryType*> METCorrectionAlgorithm::type2Binning_
private