CMS 3D CMS Logo

SignAlgoResolutions.cc
Go to the documentation of this file.
2 // -*- C++ -*-
3 //
4 // Package: METAlgorithms
5 // Class: SignAlgoResolutions
6 //
14 //
15 // Original Author: Kyle Story, Freya Blekman (Cornell University)
16 // Created: Fri Apr 18 11:58:33 CEST 2008
17 //
18 //
24 
25 #include <cmath>
26 
27 #include <cstdlib>
28 #include <iostream>
29 #include <sstream>
30 #include <string>
31 #include <sys/stat.h>
32 
34  : functionmap_(), ptResol_(nullptr), phiResol_(nullptr) {
35  addResolutions(iConfig);
36 }
37 
39  const resolutionFunc &func,
40  const double &et,
41  const double &phi,
42  const double &eta) const {
43  // derive p from et and eta;
44  double theta = 2 * atan(exp(-eta));
45  double p = et / sin(theta); // rough assumption: take e and p equivalent...
46  return eval(type, func, et, phi, eta, p);
47 }
48 
50  const resolutionFunc &func,
51  const double &et,
52  const double &phi,
53  const double &eta,
54  const double &p) const {
55  functionPars x(4);
56  x[0] = et;
57  x[1] = phi;
58  x[2] = eta;
59  x[3] = p;
60  // std::cout << "getting function of type " << type << " " << func << " " << x[0] << " " << x[1] << " " << x[2] << " " << x[3] << std::endl;
61  return getfunc(type, func, x);
62 }
63 
65  double eta = candidate->eta();
66  double phi = candidate->phi();
67  double et = candidate->energy() * sin(candidate->theta());
68  resolutionType thetype;
70  int type = candidate->particleId();
71  switch (type) {
72  case 1:
73  thetype = PFtype1;
74  name = "PFChargedHadron";
75  break;
76  case 2:
77  thetype = PFtype2;
78  name = "PFChargedEM";
79  break;
80  case 3:
81  thetype = PFtype3;
82  name = "PFMuon";
83  break;
84  case 4:
85  thetype = PFtype4;
86  name = "PFNeutralEM";
87  break;
88  case 5:
89  thetype = PFtype5;
90  name = "PFNeutralHadron";
91  break;
92  case 6:
93  thetype = PFtype6;
94  name = "PFtype6";
95  break;
96  case 7:
97  thetype = PFtype7;
98  name = "PFtype7";
99  break;
100  default:
101  thetype = PFtype7;
102  name = "PFunknown";
103  break;
104  }
105 
106  double d_et = 0, d_phi = 0; //d_phi here is the error on phi component of the et
107  reco::TrackRef trackRef = candidate->trackRef();
108  if (!trackRef.isNull()) {
109  d_phi = et * trackRef->phiError();
110  d_et = (type == 2) ? ElectronPtResolution(candidate) : trackRef->ptError();
111  //if(type==2) std::cout << eval(thetype,ET,et,phi,eta) << " " << trackRef->ptError() << " "<< ElectronPtResolution(candidate) << std::endl;
112  } else {
113  d_et = eval(thetype, ET, et, phi, eta);
114  d_phi = eval(thetype, PHI, et, phi, eta);
115  }
116 
117  metsig::SigInputObj resultingobj(name, et, phi, d_et, d_phi);
118  return resultingobj;
119 }
120 
122  double jpt = jet->pt();
123  double jphi = jet->phi();
124  double jeta = jet->eta();
125  double jdeltapt = 999.;
126  double jdeltapphi = 999.;
127 
128  if (jpt < ptResolThreshold_ && jpt < 20.) { //use temporary fix for low pT jets
129  double feta = TMath::Abs(jeta);
130  int ieta = feta < 5. ? int(feta / 0.5) : 9; //bin size = 0.5
131  int ipt = jpt > 3. ? int(jpt - 3. / 2) : 0; //bin size =2, starting from ptmin=3GeV
132  jdeltapt = jdpt[ieta][ipt];
133  jdeltapphi = jpt * jdphi[ieta][ipt];
134  } else {
135  //use the resolution functions at |eta|=5 to avoid crash for jets with large eta.
136  if (jeta > 5)
137  jeta = 5;
138  if (jeta < -5)
139  jeta = -5;
140  double jptForEval = jpt > ptResolThreshold_ ? jpt : ptResolThreshold_;
141  jdeltapt = jpt * ptResol_->parameterEtaEval("sigma", jeta, jptForEval);
142  jdeltapphi = jpt * phiResol_->parameterEtaEval("sigma", jeta, jptForEval);
143  }
144 
145  std::string inputtype = "jet";
146  metsig::SigInputObj obj_jet(inputtype, jpt, jphi, jdeltapt, jdeltapphi);
147  //std::cout << "RESOLUTIONS JET: " << jpt << " " << jphi<< " " <<jdeltapt << " " << jdeltapphi << std::endl;
148  return obj_jet;
149 }
150 
152  using namespace std;
153 
154  // Jet Resolutions - for now load from the files. Migrate to EventSetup asap.
156 
157  ptResolThreshold_ = iConfig.getParameter<double>("ptresolthreshold");
158 
159  //get temporary low pT pfjet resolutions
160  for (int ieta = 0; ieta < 10; ieta++) {
161  jdpt[ieta] = iConfig.getParameter<std::vector<double> >(Form("jdpt%d", ieta));
162  jdphi[ieta] = iConfig.getParameter<std::vector<double> >(Form("jdphi%d", ieta));
163  }
164 
165  // for now: do this by hand - this can obviously also be done via ESSource etc.
166  functionPars etparameters(3, 0);
167  functionPars phiparameters(1, 0);
168  // set the parameters per function:
169  // ECAL, BARREL:
170  std::vector<double> ebet = iConfig.getParameter<std::vector<double> >("EB_EtResPar");
171  std::vector<double> ebphi = iConfig.getParameter<std::vector<double> >("EB_PhiResPar");
172 
173  etparameters[0] = ebet[0];
174  etparameters[1] = ebet[1];
175  etparameters[2] = ebet[2];
176  phiparameters[0] = ebphi[0];
177  addfunction(caloEB, ET, etparameters);
178  addfunction(caloEB, PHI, phiparameters);
179  // ECAL, ENDCAP:
180  std::vector<double> eeet = iConfig.getParameter<std::vector<double> >("EE_EtResPar");
181  std::vector<double> eephi = iConfig.getParameter<std::vector<double> >("EE_PhiResPar");
182 
183  etparameters[0] = eeet[0];
184  etparameters[1] = eeet[1];
185  etparameters[2] = eeet[2];
186  phiparameters[0] = eephi[0];
187  addfunction(caloEE, ET, etparameters);
188  addfunction(caloEE, PHI, phiparameters);
189  // HCAL, BARREL:
190  std::vector<double> hbet = iConfig.getParameter<std::vector<double> >("HB_EtResPar");
191  std::vector<double> hbphi = iConfig.getParameter<std::vector<double> >("HB_PhiResPar");
192 
193  etparameters[0] = hbet[0];
194  etparameters[1] = hbet[1];
195  etparameters[2] = hbet[2];
196  phiparameters[0] = hbphi[0];
197  addfunction(caloHB, ET, etparameters);
198  addfunction(caloHB, PHI, phiparameters);
199  // HCAL, ENDCAP:
200  std::vector<double> heet = iConfig.getParameter<std::vector<double> >("HE_EtResPar");
201  std::vector<double> hephi = iConfig.getParameter<std::vector<double> >("HE_PhiResPar");
202 
203  etparameters[0] = heet[0];
204  etparameters[1] = heet[1];
205  etparameters[2] = heet[2];
206  phiparameters[0] = hephi[0];
207  addfunction(caloHE, ET, etparameters);
208  addfunction(caloHE, PHI, phiparameters);
209  // HCAL, Outer
210  std::vector<double> hoet = iConfig.getParameter<std::vector<double> >("HO_EtResPar");
211  std::vector<double> hophi = iConfig.getParameter<std::vector<double> >("HO_PhiResPar");
212 
213  etparameters[0] = hoet[0];
214  etparameters[1] = hoet[1];
215  etparameters[2] = hoet[2];
216  phiparameters[0] = hophi[0];
217  addfunction(caloHO, ET, etparameters);
218  addfunction(caloHO, PHI, phiparameters);
219  // HCAL, Forward
220  std::vector<double> hfet = iConfig.getParameter<std::vector<double> >("HF_EtResPar");
221  std::vector<double> hfphi = iConfig.getParameter<std::vector<double> >("HF_PhiResPar");
222 
223  etparameters[0] = hfet[0];
224  etparameters[1] = hfet[1];
225  etparameters[2] = hfet[2];
226  phiparameters[0] = hfphi[0];
227  addfunction(caloHF, ET, etparameters);
228  addfunction(caloHF, PHI, phiparameters);
229 
230  // PF objects:
231  // type 1:
232  std::vector<double> pf1et = iConfig.getParameter<std::vector<double> >("PF_EtResType1");
233  std::vector<double> pf1phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType1");
234  etparameters[0] = pf1et[0];
235  etparameters[1] = pf1et[1];
236  etparameters[2] = pf1et[2];
237  phiparameters[0] = pf1phi[0];
238  addfunction(PFtype1, ET, etparameters);
239  addfunction(PFtype1, PHI, phiparameters);
240 
241  // PF objects:
242  // type 2:
243  std::vector<double> pf2et = iConfig.getParameter<std::vector<double> >("PF_EtResType2");
244  std::vector<double> pf2phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType2");
245  etparameters[0] = pf2et[0];
246  etparameters[1] = pf2et[1];
247  etparameters[2] = pf2et[2];
248  phiparameters[0] = pf2phi[0];
249  addfunction(PFtype2, ET, etparameters);
250  addfunction(PFtype2, PHI, phiparameters);
251 
252  // PF objects:
253  // type 3:
254  std::vector<double> pf3et = iConfig.getParameter<std::vector<double> >("PF_EtResType3");
255  std::vector<double> pf3phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType3");
256  etparameters[0] = pf3et[0];
257  etparameters[1] = pf3et[1];
258  etparameters[2] = pf3et[2];
259  phiparameters[0] = pf3phi[0];
260  addfunction(PFtype3, ET, etparameters);
261  addfunction(PFtype3, PHI, phiparameters);
262 
263  // PF objects:
264  // type 4:
265  std::vector<double> pf4et = iConfig.getParameter<std::vector<double> >("PF_EtResType4");
266  std::vector<double> pf4phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType4");
267  etparameters[0] = pf4et[0];
268  etparameters[1] = pf4et[1];
269  etparameters[2] = pf4et[2];
270  //phiparameters[0]=pf4phi[0];
271  addfunction(PFtype4, ET, etparameters);
272  addfunction(PFtype4, PHI, pf4phi); //use the same functional form for photon phi error as for pT, pass whole vector
273 
274  // PF objects:
275  // type 5:
276  std::vector<double> pf5et = iConfig.getParameter<std::vector<double> >("PF_EtResType5");
277  std::vector<double> pf5phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType5");
278  etparameters[0] = pf5et[0];
279  etparameters[1] = pf5et[1];
280  etparameters[2] = pf5et[2];
281  phiparameters[0] = pf5phi[0];
282  addfunction(PFtype5, ET, etparameters);
283  addfunction(PFtype5, PHI, pf5phi);
284 
285  // PF objects:
286  // type 6:
287  std::vector<double> pf6et = iConfig.getParameter<std::vector<double> >("PF_EtResType6");
288  std::vector<double> pf6phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType6");
289  etparameters[0] = pf6et[0];
290  etparameters[1] = pf6et[1];
291  etparameters[2] = pf6et[2];
292  phiparameters[0] = pf6phi[0];
293  addfunction(PFtype6, ET, etparameters);
294  addfunction(PFtype6, PHI, phiparameters);
295 
296  // PF objects:
297  // type 7:
298  std::vector<double> pf7et = iConfig.getParameter<std::vector<double> >("PF_EtResType7");
299  std::vector<double> pf7phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType7");
300  etparameters[0] = pf7et[0];
301  etparameters[1] = pf7et[1];
302  etparameters[2] = pf7et[2];
303  phiparameters[0] = pf7phi[0];
304  addfunction(PFtype7, ET, etparameters);
305  addfunction(PFtype7, PHI, phiparameters);
306 
307  return;
308 }
309 
312  const functionPars &parameters) {
313  functionCombo mypair(type, func);
314  functionmap_[mypair] = parameters;
315 }
316 
319  functionPars &x) const {
320  double result = 0;
321  functionCombo mypair(type, func);
322 
323  if (functionmap_.count(mypair) == 0) {
324  return result;
325  }
326 
327  functionPars values = (functionmap_.find(mypair))->second;
328  switch (func) {
329  case metsig::ET:
330  return EtFunction(x, values);
331  case metsig::PHI:
332  return PhiFunction(x, values);
333  case metsig::TRACKP:
334  return PFunction(x, values);
335  case metsig::CONSTPHI:
336  return PhiConstFunction(x, values);
337  }
338 
339  // std::cout << "returning function " << type << " " << func << " " << result << " " << x[0] << std::endl;
340 
341  return result;
342 }
343 
345  if (par.size() < 3)
346  return 0.;
347  if (x.empty())
348  return 0.;
349  double et = x[0];
350  if (et <= 0.)
351  return 0.;
352  double result = et * sqrt((par[2] * par[2]) + (par[1] * par[1] / et) + (par[0] * par[0] / (et * et)));
353  return result;
354 }
355 
357  double et = x[0];
358  if (et <= 0.) {
359  return 0.;
360  }
361 
362  //if 1 parameter is C provided, returns C*pT, if three parameters N, S, C are provided, it returns the usual resolution value, as for sigmaPt
363  if (par.size() != 1 && par.size() != 3) { //only 1 or 3 parameters supported for phi function
364  return 0.;
365  } else if (par.size() == 1) {
366  return par[0] * et;
367  } else {
368  return et * sqrt((par[2] * par[2]) + (par[1] * par[1] / et) + (par[0] * par[0] / (et * et)));
369  }
370 }
372  // not currently implemented
373  return 0;
374 }
375 
377  return par[0];
378 }
379 
381  using namespace std;
382 
383  // only reinitialize the resolutsion if the pointers are zero
384  if (ptResol_ == nullptr) {
385  string resolutionsAlgo = iConfig.getParameter<std::string>("resolutionsAlgo");
386  string resolutionsEra = iConfig.getParameter<std::string>("resolutionsEra");
387 
388  string cmssw_base(std::getenv("CMSSW_BASE"));
389  string cmssw_release_base(std::getenv("CMSSW_RELEASE_BASE"));
390  string path = cmssw_base + "/src/CondFormats/JetMETObjects/data";
391  struct stat st;
392  if (stat(path.c_str(), &st) != 0) {
393  path = cmssw_release_base + "/src/CondFormats/JetMETObjects/data";
394  }
395  if (stat(path.c_str(), &st) != 0) {
396  cerr << "ERROR: tried to set path but failed, abort." << endl;
397  }
398  const string &era(resolutionsEra);
399  const string &alg(resolutionsAlgo);
400  string ptFileName = path + "/" + era + "_PtResolution_" + alg + ".txt";
401  string phiFileName = path + "/" + era + "_PhiResolution_" + alg + ".txt";
402 
403  ptResol_ = new JetResolution(ptFileName, false);
404  phiResol_ = new JetResolution(phiFileName, false);
405  }
406 }
407 
409  double eta = c->eta();
410  double energy = c->energy();
411  double dEnergy = pfresol_->getEnergyResolutionEm(energy, eta);
412 
413  return dEnergy / cosh(eta);
414 }
reco::PFCandidate::trackRef
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:408
metsig::caloHB
Definition: SignAlgoResolutions.h:46
metsig::caloHE
Definition: SignAlgoResolutions.h:43
BeamSpotPI::parameters
parameters
Definition: BeamSpotPayloadInspectorHelper.h:30
metsig::SignAlgoResolutions::PhiFunction
double PhiFunction(const functionPars &x, const functionPars &par) const
Definition: SignAlgoResolutions.cc:356
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
metsig::TRACKP
Definition: SignAlgoResolutions.h:59
metsig::SignAlgoResolutions::PFunction
double PFunction(const functionPars &x, const functionPars &par) const
Definition: SignAlgoResolutions.cc:371
metsig::ET
Definition: SignAlgoResolutions.h:59
metsig::SignAlgoResolutions::evalPF
metsig::SigInputObj evalPF(const reco::PFCandidate *candidate) const
Definition: SignAlgoResolutions.cc:64
edm::Ref::isNull
bool isNull() const
Checks for null.
Definition: Ref.h:235
metsig::SignAlgoResolutions::PhiConstFunction
double PhiConstFunction(const functionPars &x, const functionPars &par) const
Definition: SignAlgoResolutions.cc:376
METSignificanceParams_cfi.jeta
jeta
Definition: METSignificanceParams_cfi.py:12
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
metsig::SignAlgoResolutions::getfunc
double getfunc(const resolutionType &type, const resolutionFunc &func, std::vector< double > &x) const
Definition: SignAlgoResolutions.cc:317
metsig::SignAlgoResolutions::SignAlgoResolutions
SignAlgoResolutions()
Definition: SignAlgoResolutions.h:63
runTheMatrix.cmssw_base
cmssw_base
Definition: runTheMatrix.py:323
metsig::SignAlgoResolutions::initializeJetResolutions
void initializeJetResolutions(const edm::ParameterSet &iConfig)
Definition: SignAlgoResolutions.cc:380
metsig::PFtype4
Definition: SignAlgoResolutions.h:54
metsig::SignAlgoResolutions::functionPars
std::vector< double > functionPars
Definition: SignAlgoResolutions.h:88
metsig::SignAlgoResolutions::eval
double eval(const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
Definition: SignAlgoResolutions.cc:49
edm::Ref< TrackCollection >
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
SignAlgoResolutions.h
metsig::SignAlgoResolutions::ElectronPtResolution
double ElectronPtResolution(const reco::PFCandidate *c) const
Definition: SignAlgoResolutions.cc:408
Track.h
TrackFwd.h
Abs
T Abs(T a)
Definition: MathUtil.h:49
metsig::caloEB
Definition: SignAlgoResolutions.h:42
PVValHelper::eta
Definition: PVValidationHelpers.h:70
contentValuesCheck.values
values
Definition: contentValuesCheck.py:38
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
metsig::PFtype3
Definition: SignAlgoResolutions.h:53
PFBlockElement.h
reco::LeafCandidate::theta
double theta() const final
momentum polar angle
Definition: LeafCandidate.h:150
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
jpt
Definition: JetPlusTrackCorrector.h:40
metsig::PFtype6
Definition: SignAlgoResolutions.h:56
metsig::caloHO
Definition: SignAlgoResolutions.h:44
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
JetMETHLTOfflineSource_cfi.feta
feta
Definition: JetMETHLTOfflineSource_cfi.py:30
edm::ParameterSet
Definition: ParameterSet.h:47
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
reco::LeafCandidate::eta
double eta() const final
momentum pseudorapidity
Definition: LeafCandidate.h:152
metsig::resolutionType
resolutionType
Definition: SignAlgoResolutions.h:40
metsig::PFtype1
Definition: SignAlgoResolutions.h:51
metsig::caloEE
Definition: SignAlgoResolutions.h:41
reco::PFCandidate::particleId
virtual ParticleType particleId() const
Definition: PFCandidate.h:367
type
type
Definition: SiPixelVCal_PayloadInspector.cc:39
JetResolution
Definition: JetResolution.h:11
createfilelist.int
int
Definition: createfilelist.py:10
metsig::SignAlgoResolutions::EtFunction
double EtFunction(const functionPars &x, const functionPars &par) const
Definition: SignAlgoResolutions.cc:344
metsig::SignAlgoResolutions::addResolutions
void addResolutions(const edm::ParameterSet &iConfig)
Definition: SignAlgoResolutions.cc:151
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
metsig::PFtype2
Definition: SignAlgoResolutions.h:52
TrackCollections2monitor_cff.func
func
Definition: TrackCollections2monitor_cff.py:359
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
metsig::SignAlgoResolutions::addfunction
void addfunction(const resolutionType type, const resolutionFunc func, const std::vector< double > &parameters)
Definition: SignAlgoResolutions.cc:310
metsig::SignAlgoResolutions::functionCombo
std::pair< metsig::resolutionType, metsig::resolutionFunc > functionCombo
Definition: SignAlgoResolutions.h:87
PFBlock.h
METSigParams_cfi.resolutionsEra
resolutionsEra
Definition: METSigParams_cfi.py:50
metsig::PHI
Definition: SignAlgoResolutions.h:59
std
Definition: JetResolutionObject.h:76
reco::LeafCandidate::phi
double phi() const final
momentum azimuthal angle
Definition: LeafCandidate.h:148
class-composition.getfunc
getfunc
Definition: class-composition.py:16
metsig::jet
Definition: SignAlgoResolutions.h:47
metsig::PFtype7
Definition: SignAlgoResolutions.h:57
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
reco::LeafCandidate::energy
double energy() const final
energy
Definition: LeafCandidate.h:125
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
metsig::resolutionFunc
resolutionFunc
Definition: SignAlgoResolutions.h:59
mps_fire.result
result
Definition: mps_fire.py:311
castor_dqm_sourceclient_file_cfg.path
path
Definition: castor_dqm_sourceclient_file_cfg.py:37
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:56
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
edm_modernize_messagelogger.stat
stat
Definition: edm_modernize_messagelogger.py:27
METSigParams_cfi.resolutionsAlgo
resolutionsAlgo
Definition: METSigParams_cfi.py:51
metsig::PFtype5
Definition: SignAlgoResolutions.h:55
EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.cerr
cerr
Definition: EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.py:8
metsig::SigInputObj
Definition: SigInputObj.h:29
metsig::caloHF
Definition: SignAlgoResolutions.h:45
metsig::SignAlgoResolutions::evalPFJet
metsig::SigInputObj evalPFJet(const reco::Jet *jet) const
Definition: SignAlgoResolutions.cc:121
metsig::CONSTPHI
Definition: SignAlgoResolutions.h:59