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 &, const edm::EventSetup &)
 
 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

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

Definition at line 43 of file METCorrectionAlgorithm.h.

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

Definition at line 55 of file METCorrectionAlgorithm.h.

Constructor & Destructor Documentation

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

Definition at line 9 of file METCorrectionAlgorithm.cc.

References applyType0Corrections_, applyType1Corrections_, applyType2Corrections_, chsSumTokens_, data, Exception, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), HiClusterCompatibility_cfi::inputTag, patCaloMETCorrections_cff::srcType1Corrections, patCaloMETCorrections_cff::srcUnclEnergySums, AlCaHLTBitMon_QueryRunRegistry::string, type0Cuncl_, type0Rsoft_, type1Tokens_, type2Binning_, patCaloMETCorrections_cff::type2CorrFormula, and patCaloMETCorrections_cff::type2CorrParameter.

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

Definition at line 65 of file METCorrectionAlgorithm.cc.

References type2Binning_.

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

Member Function Documentation

CorrMETData METCorrectionAlgorithm::compMETCorrection ( edm::Event evt,
const edm::EventSetup es 
)

Definition at line 73 of file METCorrectionAlgorithm.cc.

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

74 {
75  CorrMETData metCorr;
76  metCorr.mex = 0.;
77  metCorr.mey = 0.;
78  metCorr.sumet = 0.;
79 
80  if ( applyType0Corrections_ ) {
81 //--- sum all Type 0 MET correction terms
83  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = chsSumTokens_.begin(); corrToken != chsSumTokens_.end(); ++corrToken)
84  {
85  evt.getByToken(*corrToken, chsSum);
86 
87  metCorr.mex += type0Cuncl_*(1 - type0Rsoft_)*chsSum->mex;
88  metCorr.mey += type0Cuncl_*(1 - type0Rsoft_)*chsSum->mey;
89  metCorr.sumet += type0Cuncl_*(1 - type0Rsoft_)*chsSum->sumet;
90  }
91  }
92 
93  if ( applyType1Corrections_ ) {
94 //--- sum all Type 1 MET correction terms
95  edm::Handle<CorrMETData> type1Correction;
96  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = type1Tokens_.begin(); corrToken != type1Tokens_.end(); ++corrToken)
97  {
98  evt.getByToken(*corrToken, type1Correction);
99 
100  metCorr.mex += type1Correction->mex;
101  metCorr.mey += type1Correction->mey;
102  metCorr.sumet += type1Correction->sumet;
103  }
104  }
105 
107  {
108 //--- compute momentum sum of all "unclustered energy" in the event
109 //
110 // NOTE: calibration factors/formulas for Type 2 MET correction may depend on eta
111 // (like the jet energy correction factors do)
112 //
113 
114  for ( std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
115  type2BinningEntry != type2Binning_.end(); ++type2BinningEntry )
116  {
117  CorrMETData unclEnergySum;
118  edm::Handle<CorrMETData> unclEnergySummand;
119  for (std::vector<edm::EDGetTokenT<CorrMETData> >::const_iterator corrToken = (*type2BinningEntry)->corrTokens_.begin(); corrToken != (*type2BinningEntry)->corrTokens_.end(); ++corrToken)
120  {
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:460
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< edm::EDGetTokenT< CorrMETData > > chsSumTokens_
double sumet
Definition: CorrMETData.h:20
std::vector< type2BinningEntryType * > type2Binning_
std::vector< edm::EDGetTokenT< CorrMETData > > type1Tokens_
a MET correction term
Definition: CorrMETData.h:14
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17

Member Data Documentation

bool METCorrectionAlgorithm::applyType0Corrections_
private

Definition at line 48 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

bool METCorrectionAlgorithm::applyType1Corrections_
private

Definition at line 49 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

bool METCorrectionAlgorithm::applyType2Corrections_
private

Definition at line 50 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

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

Definition at line 45 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

double METCorrectionAlgorithm::type0Cuncl_
private

Definition at line 53 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

double METCorrectionAlgorithm::type0Rsoft_
private

Definition at line 52 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

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

Definition at line 46 of file METCorrectionAlgorithm.h.

Referenced by compMETCorrection(), and METCorrectionAlgorithm().

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