CMS 3D CMS Logo

EGEtScaleSysModifier.cc
Go to the documentation of this file.
5 
6 #include <vdt/vdtMath.h>
7 
8 //in the legacy re-reco we came across an Et scale issue, specifically an inflection at 45 GeV
9 //it is problematic to address in the normal scale and smearing code therefore
10 //this patch modifies the e/gamma object to "patch in" this additional systematic
11 //Questions:
12 //1 ) Why dont we add this to the overall scale systematic?
13 //Well we could but the whole point of this is to get a proper template which can be evolved
14 //correctly. The issue is not that the current systematic doesnt cover this systematic on a
15 //per bin basis, it does, the problem is the sign flips at 45 GeV so its fine if you use
16 //only eles/phos > 45 or <45 GeV but the template cant simultaneously cover the entire spectrum
17 //and you'll get a nasty constrained fit.
18 //And if you care about these sort of things (and you probably should), you shouldnt be using
19 //the overall systematic anyways but its seperate parts
20 //
21 //2 ) Why dont you do it inside the standard class
22 //We could but we're rather hoping to solve this issue more cleanly in the future
23 //but we would hold up the reminiAOD too long to do so so this is just to get something
24 //which should be fine for 90% of analyses in the miniAOD. So this is a temporary fix
25 //which we hope to rip out soon (if you're reading this in 2024 because we're still doing it
26 //this way, all I can say is sorry, we thought it would go away soon! )
27 
28 
30 public:
32  ~EGEtScaleSysModifier() override{}
33 
34  void setEvent(const edm::Event&) final;
35  void setEventContent(const edm::EventSetup&) final;
37 
38  void modifyObject(pat::Electron& ele) const final;
39  void modifyObject(pat::Photon& pho) const final;
40 
41 private:
42  std::pair<float,float> calCombinedMom(reco::GsfElectron& ele,const float scale,
43  const float smear)const;
44  void setEcalEnergy(reco::GsfElectron& ele,const float scale,const float smear)const;
45 
47  public:
49  virtual ~UncertFuncBase(){}
50  virtual float val(const float et)const=0;
51  };
52 
53  //defines two uncertaintes, one Et<X and one Et>Y
54  //for X<Et<Y, it linearly extrapolates betwen the two values
55  class UncertFuncV1 : public UncertFuncBase {
56  public:
58  lowEt_(conf.getParameter<double>("lowEt")),highEt_(conf.getParameter<double>("highEt")),
59  lowEtUncert_(conf.getParameter<double>("lowEtUncert")),
60  highEtUncert_(conf.getParameter<double>("highEtUncert")),
61  dEt_(highEt_-lowEt_),
62  dUncert_(highEtUncert_ - lowEtUncert_)
63  {
64  if(highEt_<=lowEt_) throw cms::Exception("ConfigError") <<" highEt "<<highEt_<<" is not higher than lowEt "<<lowEt_;
65  }
66  ~UncertFuncV1() override{}
67 
68  float val(const float et)const override{
69  if(et<=lowEt_) return lowEtUncert_;
70  else if(et>=highEt_) return highEtUncert_;
71  else{
72  return (et-lowEt_)*dUncert_/dEt_+lowEtUncert_;
73  }
74  }
75  private:
76  float lowEt_;
77  float highEt_;
78  float lowEtUncert_;
80  float dEt_;
81  float dUncert_;
82  };
83 
85  std::unique_ptr<UncertFuncBase> uncertFunc_;
86 };
87 
90  epCombTool_(conf.getParameter<edm::ParameterSet>("epCombConfig"))
91 {
92  const edm::ParameterSet funcPSet = conf.getParameter<edm::ParameterSet>("uncertFunc");
93  const std::string& funcName = funcPSet.getParameter<std::string>("name");
94  if(funcName == "UncertFuncV1"){
95  uncertFunc_ = std::make_unique<UncertFuncV1>(funcPSet);
96  }else{
97  throw cms::Exception("ConfigError") <<"Error constructing EGEtScaleSysModifier, function name "<<funcName<<" not valid";
98  }
99 
100 }
101 
103 {
104 
105 }
106 
108 {
109 
110 }
111 
113 {
115 }
116 
117 
119 {
120  auto getVal = [](const pat::Electron& ele,EGEnergySysIndex::Index valIndex){
121  return ele.userFloat(EGEnergySysIndex::name(valIndex));
122  };
123  //so ele.energy() may be either pre or post corrections, we have no idea
124  //so we explicity access the pre and post correction ecal energies
125  //we need the pre corrected to properly do the e/p combination
126  //we need the post corrected to get et uncertainty
127  const float ecalEnergyPostCorr = getVal(ele,EGEnergySysIndex::kEcalPostCorr);
128  const float ecalEnergyPreCorr = getVal(ele,EGEnergySysIndex::kEcalPreCorr);
129  const float ecalEnergyErrPreCorr = getVal(ele,EGEnergySysIndex::kEcalErrPreCorr);
130 
131  //the et cut is in terms of ecal et using the track angle and post corr ecal energy
132  const float etUncert = uncertFunc_->val(ele.et()/ele.energy()*ecalEnergyPostCorr);
133  const float smear = getVal(ele,EGEnergySysIndex::kSmearValue);
134  const float corr = getVal(ele,EGEnergySysIndex::kScaleValue);
135 
136  //get the values we have to reset back to
137  const float oldEcalEnergy = ele.ecalEnergy();
138  const float oldEcalEnergyErr = ele.ecalEnergyError();
139  const auto oldP4 = ele.p4();
140  const float oldP4Err = ele.p4Error(reco::GsfElectron::P4_COMBINATION);
141  const float oldTrkMomErr = ele.trackMomentumError();
142 
143  ele.setCorrectedEcalEnergy(ecalEnergyPreCorr);
144  ele.setCorrectedEcalEnergyError(ecalEnergyErrPreCorr);
145 
146  const float energyEtUncertUp = calCombinedMom(ele,corr+etUncert,smear).first;
147  const float energyEtUncertDn = calCombinedMom(ele,corr-etUncert,smear).first;
148 
149  //reset it back to how it was
150  ele.setCorrectedEcalEnergy(oldEcalEnergy);
151  ele.setCorrectedEcalEnergyError(oldEcalEnergyErr);
152  ele.correctMomentum(oldP4,oldTrkMomErr,oldP4Err);
153 
154  ele.addUserFloat("energyScaleEtUp",energyEtUncertUp);
155  ele.addUserFloat("energyScaleEtDown",energyEtUncertDn);
156 
157 }
158 
160 {
161  auto getVal = [](const pat::Photon& pho,EGEnergySysIndex::Index valIndex){
162  return pho.userFloat(EGEnergySysIndex::name(valIndex));
163  };
164  //so pho.energy() may be either pre or post corrections, we have no idea
165  //so we explicity access the pre and post correction ecal energies
166  //post corr for the et value for the systematic, pre corr to apply them
167  const float ecalEnergyPostCorr = getVal(pho,EGEnergySysIndex::kEcalPostCorr);
168  const float ecalEnergyPreCorr = getVal(pho,EGEnergySysIndex::kEcalPreCorr);
169 
170  //the et cut is in terms of post corr ecal energy
171  const float etUncert = uncertFunc_->val(pho.et()/pho.energy()*ecalEnergyPostCorr);
172  const float corr = getVal(pho,EGEnergySysIndex::kScaleValue);
173 
174  const float energyEtUncertUp = ecalEnergyPreCorr*(corr+etUncert);
175  const float energyEtUncertDn = ecalEnergyPreCorr*(corr-etUncert);
176 
177  pho.addUserFloat("energyScaleEtUp",energyEtUncertUp);
178  pho.addUserFloat("energyScaleEtDown",energyEtUncertDn);
179 
180 }
181 
183  const float scale,
184  const float smear)const
185 {
186  const float oldEcalEnergy = ele.ecalEnergy();
187  const float oldEcalEnergyErr = ele.ecalEnergyError();
188  const auto oldP4 = ele.p4();
189  const float oldP4Err = ele.p4Error(reco::GsfElectron::P4_COMBINATION);
190  const float oldTrkMomErr = ele.trackMomentumError();
191 
192  setEcalEnergy(ele,scale,smear);
193  const auto& combinedMomentum = epCombTool_.combine(ele);
194  ele.setCorrectedEcalEnergy(oldEcalEnergy);
195  ele.setCorrectedEcalEnergyError(oldEcalEnergyErr);
196  ele.correctMomentum(oldP4,oldTrkMomErr,oldP4Err);
197 
198 
199  return combinedMomentum;
200 }
201 
203  const float scale,
204  const float smear)const
205 {
206  const float oldEcalEnergy = ele.ecalEnergy();
207  const float oldEcalEnergyErr = ele.ecalEnergyError();
208  ele.setCorrectedEcalEnergy( oldEcalEnergy*scale );
209  ele.setCorrectedEcalEnergyError(std::hypot( oldEcalEnergyErr*scale, oldEcalEnergy*smear*scale ) );
210 }
211 
212 
213 
216  "EGEtScaleSysModifier");
EGEtScaleSysModifier(const edm::ParameterSet &conf)
T getParameter(std::string const &) const
Analysis-level Photon class.
Definition: Photon.h:47
virtual float val(const float et) const =0
float trackMomentumError() const
Definition: GsfElectron.h:825
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:225
std::unique_ptr< UncertFuncBase > uncertFunc_
void modifyObject(pat::Electron &ele) const final
static const std::string & name(size_t index)
void correctMomentum(const LorentzVector &p4, float trackMomentumError, float p4Error)
Definition: GsfElectron.h:847
EpCombinationTool epCombTool_
void combine(SimpleElectron &mySimpleElectron) const
float p4Error(P4Kind kind) const
Definition: GsfElectron.cc:237
void addUserFloat(const std::string &label, float data, const bool overwrite=false)
Set user-defined float.
Definition: PATObject.h:813
float userFloat(const std::string &key) const
Definition: PATObject.h:791
std::pair< float, float > calCombinedMom(reco::GsfElectron &ele, const float scale, const float smear) const
int iEvent
Definition: GenABIO.cc:230
void setCorrectedEcalEnergyError(float newEnergyError)
Definition: GsfElectron.cc:179
void setConsumes(edm::ConsumesCollector &) final
double et() const final
transverse energy
double energy() const final
energy
UncertFuncV1(const edm::ParameterSet &conf)
void setEventContent(const edm::EventSetup &iSetup)
float ecalEnergyError() const
Definition: GsfElectron.h:838
JetCorrectorParameters corr
Definition: classes.h:5
void setEvent(const edm::Event &) final
Analysis-level electron class.
Definition: Electron.h:52
void setEventContent(const edm::EventSetup &) final
void setEcalEnergy(reco::GsfElectron &ele, const float scale, const float smear) const
et
define resolution functions of each parameter
float ecalEnergy() const
Definition: GsfElectron.h:837
void setCorrectedEcalEnergy(float newEnergy)
Definition: GsfElectron.cc:182
HLT enums.
#define DEFINE_EDM_PLUGIN(factory, type, name)
float val(const float et) const override