CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
METCorrectionAlgorithm.cc
Go to the documentation of this file.
2 
4 
5 #include <TString.h>
6 
7 #include <string>
8 
10 {
11  applyType1Corrections_ = cfg.getParameter<bool>("applyType1Corrections");
12  if ( applyType1Corrections_ ) {
13  srcType1Corrections_ = cfg.getParameter<vInputTag>("srcType1Corrections");
14  }
15 
16  applyType2Corrections_ = cfg.getParameter<bool>("applyType2Corrections");
17  if ( applyType2Corrections_ ) {
18  vInputTag srcUnclEnergySums = cfg.getParameter<vInputTag>("srcUnclEnergySums");
19 
20  if ( cfg.exists("type2Binning") ) {
21  typedef std::vector<edm::ParameterSet> vParameterSet;
22  vParameterSet cfgType2Binning = cfg.getParameter<vParameterSet>("type2Binning");
23  for ( vParameterSet::const_iterator cfgType2BinningEntry = cfgType2Binning.begin();
24  cfgType2BinningEntry != cfgType2Binning.end(); ++cfgType2BinningEntry ) {
25  type2Binning_.push_back(new type2BinningEntryType(*cfgType2BinningEntry, srcUnclEnergySums));
26  }
27  } else {
28  std::string type2CorrFormula = cfg.getParameter<std::string>("type2CorrFormula").data();
29  edm::ParameterSet type2CorrParameter = cfg.getParameter<edm::ParameterSet>("type2CorrParameter");
30  type2Binning_.push_back(new type2BinningEntryType(type2CorrFormula, type2CorrParameter, srcUnclEnergySums));
31  }
32  }
33 
34  applyType0Corrections_ = cfg.exists("applyType0Corrections") ? cfg.getParameter<bool>("applyType0Corrections") : false;
35  if ( applyType0Corrections_ ) {
36  srcCHSSums_ = cfg.getParameter<vInputTag>("srcCHSSums");
37  type0Rsoft_ = cfg.getParameter<double>("type0Rsoft");
38  type0Cuncl_ = 1.0;
40  if (cfg.exists("type2Binning")) throw cms::Exception("Invalid Arg") << "Currently, applyType0Corrections and type2Binning cannot be used together!";
41  std::string type2CorrFormula = cfg.getParameter<std::string>("type2CorrFormula").data();
42  if (!(type2CorrFormula == "A")) throw cms::Exception("Invalid Arg") << "type2CorrFormula must be \"A\" if applyType0Corrections!";
43  edm::ParameterSet type2CorrParameter = cfg.getParameter<edm::ParameterSet>("type2CorrParameter");
44  type0Cuncl_ = type2CorrParameter.getParameter<double>("A");
45  }
46  }
47 }
48 
50 {
51  for ( std::vector<type2BinningEntryType*>::const_iterator it = type2Binning_.begin();
52  it != type2Binning_.end(); ++it ) {
53  delete (*it);
54  }
55 }
56 
58 {
59  CorrMETData metCorr;
60  metCorr.mex = 0.;
61  metCorr.mey = 0.;
62  metCorr.sumet = 0.;
63 
64  if ( applyType0Corrections_ ) {
65 //--- sum all Type 0 MET correction terms
66  for ( vInputTag::const_iterator srcCHSSum = srcCHSSums_.begin();
67  srcCHSSum != srcCHSSums_.end(); ++srcCHSSum ) {
69  evt.getByLabel(*srcCHSSum, chsSum);
70 
71  metCorr.mex += type0Cuncl_*(1 - type0Rsoft_)*chsSum->mex;
72  metCorr.mey += type0Cuncl_*(1 - type0Rsoft_)*chsSum->mey;
73  metCorr.sumet += type0Cuncl_*(1 - type0Rsoft_)*chsSum->sumet;
74  }
75  }
76 
77  if ( applyType1Corrections_ ) {
78 //--- sum all Type 1 MET correction terms
79  for ( vInputTag::const_iterator srcType1Correction = srcType1Corrections_.begin();
80  srcType1Correction != srcType1Corrections_.end(); ++srcType1Correction ) {
81  edm::Handle<CorrMETData> type1Correction;
82  evt.getByLabel(*srcType1Correction, type1Correction);
83 
84  metCorr.mex += type1Correction->mex;
85  metCorr.mey += type1Correction->mey;
86  metCorr.sumet += type1Correction->sumet;
87  }
88  }
89 
90  if ( applyType2Corrections_ ) {
91 //--- compute momentum sum of all "unclustered energy" in the event
92 //
93 // NOTE: calibration factors/formulas for Type 2 MET correction may depend on eta
94 // (like the jet energy correction factors do)
95 //
96  for ( std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
97  type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) {
98  CorrMETData unclEnergySum;
99  for ( vInputTag::const_iterator srcUnclEnergySum = (*type2BinningEntry)->srcUnclEnergySums_.begin();
100  srcUnclEnergySum != (*type2BinningEntry)->srcUnclEnergySums_.end(); ++srcUnclEnergySum ) {
101  edm::Handle<CorrMETData> unclEnergySummand;
102  evt.getByLabel(*srcUnclEnergySum, unclEnergySummand);
103 
104  unclEnergySum.mex += unclEnergySummand->mex;
105  unclEnergySum.mey += unclEnergySummand->mey;
106  unclEnergySum.sumet += unclEnergySummand->sumet;
107  }
108 
109 //--- calibrate "unclustered energy"
110  double unclEnergySumPt = sqrt(unclEnergySum.mex*unclEnergySum.mex + unclEnergySum.mey*unclEnergySum.mey);
111  double unclEnergyScaleFactor = (*type2BinningEntry)->binCorrFormula_->Eval(unclEnergySumPt);
112 
113 //--- MET balances momentum of reconstructed particles,
114 // hence correction to "unclustered energy" and corresponding Type 2 MET correction are of opposite sign
115  metCorr.mex -= (unclEnergyScaleFactor - 1.)*unclEnergySum.mex;
116  metCorr.mey -= (unclEnergyScaleFactor - 1.)*unclEnergySum.mey;
117  metCorr.sumet += (unclEnergyScaleFactor - 1.)*unclEnergySum.sumet;
118  }
119  }
120 
121  return metCorr;
122 }
T getParameter(std::string const &) const
METCorrectionAlgorithm(const edm::ParameterSet &)
CorrMETData compMETCorrection(edm::Event &, const edm::EventSetup &)
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< edm::InputTag > vInputTag
T sqrt(T t)
Definition: SSEVec.h:46
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double sumet
Definition: CorrMETData.h:20
std::vector< type2BinningEntryType * > type2Binning_
a MET correction term
Definition: CorrMETData.h:14
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double mey
Definition: CorrMETData.h:18
double mex
Definition: CorrMETData.h:17