CMS 3D CMS Logo

MET.cc
Go to the documentation of this file.
1 //
2 //
3 
5 
6 
7 using namespace pat;
8 
9 
12  initCorMap();
13 }
14 
15 
17 MET::MET(const reco::MET & aMET) : PATObject<reco::MET>(aMET) {
18  const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(&aMET);
19  if (calo != nullptr) caloMET_.push_back(calo->getSpecific());
20  const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(&aMET);
21  if (pf != nullptr) pfMET_.push_back(pf->getSpecific());
22  const pat::MET * pm = dynamic_cast<const pat::MET *>(&aMET);
23  if (pm != nullptr) this->operator=(*pm);
24 
25  metSig_ =0.;
26  initCorMap();
27 }
28 
29 
31 MET::MET(const edm::RefToBase<reco::MET> & aMETRef) : PATObject<reco::MET>(aMETRef) {
32  const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(aMETRef.get());
33  if (calo != nullptr) caloMET_.push_back(calo->getSpecific());
34  const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(aMETRef.get());
35  if (pf != nullptr) pfMET_.push_back(pf->getSpecific());
36  const pat::MET * pm = dynamic_cast<const pat::MET *>(aMETRef.get());
37  if (pm != nullptr) this->operator=(*pm);
38 
39  metSig_ =0.;
40  initCorMap();
41 }
42 
44 MET::MET(const edm::Ptr<reco::MET> & aMETRef) : PATObject<reco::MET>(aMETRef) {
45  const reco::CaloMET * calo = dynamic_cast<const reco::CaloMET *>(aMETRef.get());
46  if (calo != nullptr) caloMET_.push_back(calo->getSpecific());
47  const reco::PFMET * pf = dynamic_cast<const reco::PFMET *>(aMETRef.get());
48  if (pf != nullptr) pfMET_.push_back(pf->getSpecific());
49  const pat::MET * pm = dynamic_cast<const pat::MET *>(aMETRef.get());
50  if (pm != nullptr) this->operator=(*pm);
51 
52  metSig_ =0.;
53  initCorMap();
54 }
55 
57 MET::MET(MET const& iOther):
58 PATObject<reco::MET>(iOther),
59 genMET_(iOther.genMET_),
60 caloMET_(iOther.caloMET_),
61 pfMET_(iOther.pfMET_),
62 metSig_(iOther.metSig_),
63 uncertaintiesRaw_(iOther.uncertaintiesRaw_), //74X reading compatibility
64 uncertaintiesType1_(iOther.uncertaintiesType1_), //74X compatibility
65 uncertaintiesType1p2_(iOther.uncertaintiesType1p2_), //74X compatibility
69 
70  initCorMap();
71 }
72 
74 // old uncertainties discarded on purpose to avoid confusion
75 MET::MET(const reco::MET & corMET, const MET& srcMET ):
76 PATObject<reco::MET>(corMET),
77 genMET_(srcMET.genMET_),
78 caloMET_(srcMET.caloMET_),
79 pfMET_(srcMET.pfMET_),
80 metSig_(srcMET.metSig_),
82 
84 
85  initCorMap();
86 }
87 
90 
91 }
92 
93 MET& MET::operator=(MET const& iOther) {
95  genMET_ = iOther.genMET_;
96  caloMET_ =iOther.caloMET_;
97  pfMET_ =iOther.pfMET_;
98  uncertaintiesRaw_ = iOther.uncertaintiesRaw_; //74X compatibility
102  corrections_ = iOther.corrections_;
103  metSig_ = iOther.metSig_;
105 
106  return *this;
107 }
108 
110 const reco::GenMET * MET::genMET() const {
111  return (!genMET_.empty() ? &genMET_.front() : nullptr );
112 }
113 
115 void MET::setGenMET(const reco::GenMET & gm) {
116  genMET_.clear();
117  genMET_.push_back(gm);
118 }
119 
120 
121 //Method to set the MET significance
122 void MET::setMETSignificance(const double& metSig) {
123  metSig_ = metSig;
124 }
125 
126 double MET::metSignificance() const {
127  return metSig_;
128 }
129 
130 void
132 
133  std::vector<MET::METCorrectionType> tmpRaw;
134  std::vector<MET::METCorrectionType> tmpType1;
135  std::vector<MET::METCorrectionType> tmpType01;
136  std::vector<MET::METCorrectionType> tmpTypeXY;
137  std::vector<MET::METCorrectionType> tmpType1XY;
138  std::vector<MET::METCorrectionType> tmpType01XY;
139  std::vector<MET::METCorrectionType> tmpType1Smear;
140  std::vector<MET::METCorrectionType> tmpType01Smear;
141  std::vector<MET::METCorrectionType> tmpType1SmearXY;
142  std::vector<MET::METCorrectionType> tmpType01SmearXY;
143 
144  tmpRaw.push_back(MET::None);
145 
146  tmpType1.push_back(MET::T1);
147  tmpType01.push_back(MET::T1);
148  tmpType1XY.push_back(MET::T1);
149  tmpType01XY.push_back(MET::T1);
150  tmpType1Smear.push_back(MET::T1);
151  tmpType01Smear.push_back(MET::T1);
152  tmpType1SmearXY.push_back(MET::T1);
153  tmpType01SmearXY.push_back(MET::T1);
154 
155  tmpType01.push_back(MET::T0);
156  tmpType01XY.push_back(MET::T0);
157  tmpType01Smear.push_back(MET::T0);
158  tmpType01SmearXY.push_back(MET::T0);
159 
160  tmpType1Smear.push_back(MET::Smear);
161  tmpType01Smear.push_back(MET::Smear);
162  tmpType1SmearXY.push_back(MET::Smear);
163  tmpType01SmearXY.push_back(MET::Smear);
164 
165  tmpTypeXY.push_back(MET::TXYForRaw);
166  tmpType1XY.push_back(MET::TXY);
167  tmpType01XY.push_back(MET::TXYForT01);
168  tmpType1SmearXY.push_back(MET::TXYForT1Smear);
169  tmpType01SmearXY.push_back(MET::TXYForT01Smear);
170 
171  corMap_[MET::Raw] = tmpRaw;
172  corMap_[MET::Type1] = tmpType1;
173  corMap_[MET::Type01] = tmpType01;
174  corMap_[MET::TypeXY] = tmpTypeXY;
175  corMap_[MET::Type1XY] = tmpType1XY;
176  corMap_[MET::Type01XY] = tmpType01XY;
177  corMap_[MET::Type1Smear] = tmpType1Smear;
178  corMap_[MET::Type01Smear] = tmpType01Smear;
179  corMap_[MET::Type1SmearXY] = tmpType1SmearXY;
180  corMap_[MET::Type01SmearXY] = tmpType01SmearXY;
181 
182  //specific calo case
183  std::vector<MET::METCorrectionType> tmpRawCalo;
184  tmpRawCalo.push_back(MET::Calo);
185  corMap_[MET::RawCalo] = tmpRawCalo;
186 
187  //specific chs case
188  std::vector<MET::METCorrectionType> tmpRawChs;
189  tmpRawChs.push_back(MET::Chs);
190  corMap_[MET::RawChs] = tmpRawChs;
191 
192  //specific trk case
193  std::vector<MET::METCorrectionType> tmpRawTrk;
194  tmpRawTrk.push_back(MET::Trk);
195  corMap_[MET::RawTrk] = tmpRawTrk;
196 }
197 
200 
201  //find corrections shifts =============================
202  std::map<MET::METCorrectionLevel, std::vector<MET::METCorrectionType> >::const_iterator itCor_ = corMap_.find(cor);
203  if(itCor_==corMap_.end() ) throw cms::Exception("Unsupported", "Specified MET correction scheme does not exist");
204 
205  bool isSmeared=false;
206  MET::PackedMETUncertainty totShift;
207  unsigned int scor=itCor_->second.size();
208  for(unsigned int i=0; i<scor;i++) {
209  totShift.add( corrections_[ itCor_->second[i] ].dpx(),
210  corrections_[ itCor_->second[i] ].dpy(),
211  corrections_[ itCor_->second[i] ].dsumEt() );
212 
213  if(itCor_->first>=MET::Type1Smear)
214  isSmeared=true;
215  }
216 
217 
218 
219  //find uncertainty shift =============================
220 
221  if (uncertainties_.empty())
222  return totShift;
223 
224  if(shift>=MET::METUncertaintySize) throw cms::Exception("Unsupported", "MET uncertainty does not exist");
225  if(isSmeared && shift<=MET::JetResDown) shift = (MET::METUncertainty)(MET::METUncertaintySize+shift+1);
226 
227  totShift.add( uncertainties_[ shift ].dpx(),
228  uncertainties_[ shift ].dpy(),
229  uncertainties_[ shift ].dsumEt() );
230 
231  return totShift;
232 }
233 
234 
236 
237  Vector2 vo;
238 
239  //backward compatibility with 74X samples -> the only one
240  // with uncertaintiesType1_/uncertaintiesRaw_ not empty
241  //will be removed once 74X is not used anymore
242  if(!uncertaintiesType1_.empty() || !uncertaintiesRaw_.empty()) {
243  if(cor!=MET::METCorrectionLevel::RawCalo) {
244  vo = shiftedP2_74x(shift, cor);
245  } else {
247  vo = ret;
248  }
249  }
250  else {
251  const MET::PackedMETUncertainty& v = findMETTotalShift(cor,shift);
252  Vector2 ret{ (px() + v.dpx()), (py() + v.dpy()) };
253  //return ret;
254  vo = ret;
255  }
256  return vo;
257 }
259 
260  Vector vo;
261 
262  //backward compatibility with 74X samples -> the only one
263  // with uncertaintiesType1_/uncertaintiesRaw_ not empty
264  //will be removed once 74X is not used anymore
265  if(!uncertaintiesType1_.empty() || !uncertaintiesRaw_.empty()) {
266  if(cor!=MET::METCorrectionLevel::RawCalo) {
267  vo = shiftedP3_74x(shift, cor);
268  } else {
270  vo = tmp;
271  }
272  }
273  else {
274  const MET::PackedMETUncertainty& v = findMETTotalShift(cor,shift);
275  //return Vector(px() + v.dpx(), py() + v.dpy(), 0);
276  Vector tmp(px() + v.dpx(), py() + v.dpy(), 0);
277  vo = tmp;
278  }
279  return vo;
280 }
282 
283  LorentzVector vo;
284 
285  //backward compatibility with 74X samples -> the only one
286  // with uncertaintiesType1_/uncertaintiesRaw_ not empty
287  //will be removed once 74X is not used anymore
288  if(!uncertaintiesType1_.empty() || !uncertaintiesRaw_.empty()) {
289  if(cor!=MET::METCorrectionLevel::RawCalo) {
290  vo = shiftedP4_74x(shift, cor);
291  } else {
292  double x = caloPackedMet_.dpx(), y = caloPackedMet_.dpy();
293  LorentzVector tmp(x, y, 0, std::hypot(x,y));
294  vo = tmp;
295  }
296  }
297  else {
298  const MET::PackedMETUncertainty& v = findMETTotalShift(cor,shift);
299  double x = px() + v.dpx(), y = py() + v.dpy();
300  //return LorentzVector(x, y, 0, std::hypot(x,y));
301  LorentzVector tmp(x, y, 0, std::hypot(x,y));
302  vo = tmp;
303  }
304  return vo;
305 }
307 
308  double sumEto;
309 
310  //backward compatibility with 74X samples -> the only one
311  // with uncertaintiesType1_/uncertaintiesRaw_ not empty
312  //will be removed once 74X is not used anymore
313  if(!uncertaintiesType1_.empty() || !uncertaintiesRaw_.empty()) {
314  if(cor!=MET::METCorrectionLevel::RawCalo) {
315  sumEto = shiftedSumEt_74x(shift, cor);
316  } else {
317  sumEto = caloPackedMet_.dsumEt();
318  }
319  }
320  else {
321  const MET::PackedMETUncertainty& v = findMETTotalShift(cor,shift);
322  //return sumEt() + v.dsumEt();
323  sumEto = sumEt() + v.dsumEt();
324  }
325  return sumEto;
326 }
327 
329  return shiftedP2(MET::NoShift, cor );
330 }
332  return shiftedP3(MET::NoShift, cor );
333 }
335  return shiftedP4(MET::NoShift, cor );
336 }
338  return shiftedSumEt(MET::NoShift, cor );
339 }
340 
342  return shiftedP2(MET::NoShift, MET::Raw );
343 }
345  return shiftedP3(MET::NoShift, MET::Raw );
346 }
348  return shiftedP4(MET::NoShift, MET::Raw );
349 }
350 double MET::uncorSumEt() const {
352 }
353 
354 
355 void MET::setUncShift(double px, double py, double sumEt, METUncertainty shift, bool isSmeared) {
356  if (uncertainties_.empty()) {
357  uncertainties_.resize(METUncertainty::METFullUncertaintySize);
358  }
359 
360  if(isSmeared && shift<=MET::JetResDown) {
361  //changing reference to only get the uncertainty shift and not the smeared one
362  // which is performed independently
363  shift = (MET::METUncertainty)(METUncertainty::METUncertaintySize+shift+1);
364  const PackedMETUncertainty& ref = uncertainties_[METUncertainty::NoShift];
365  uncertainties_[shift].set(px + ref.dpx() - this->px(), py + ref.dpy() - this->py(), sumEt + ref.dsumEt() - this->sumEt() );
366  }
367  else
368  uncertainties_[shift].set(px - this->px(), py - this->py(), sumEt - this->sumEt());
369 
370 }
371 
372 void MET::setCorShift(double px, double py, double sumEt, MET::METCorrectionType level) {
373  if (corrections_.empty()) {
374  corrections_.resize(MET::METCorrectionType::METCorrectionTypeSize);
375  }
376 
377  corrections_[level].set(px - this->px(), py - this->py(), sumEt - this->sumEt());
378 }
379 
380 
382  return shiftedP2(MET::METUncertainty::NoShift, MET::METCorrectionLevel::RawCalo );
383 }
384 
385 double MET::caloMETPt() const {
386  return caloMETP2().pt();
387 }
388 
389 double MET::caloMETPhi() const {
390  return caloMETP2().phi();
391 }
392 
393 double MET::caloMETSumEt() const {
395 }
396 
397 // functions to access to 74X samples ========================================================
399  if (level != Type1 && level != Raw) throw cms::Exception("Unsupported", "MET uncertainties only supported for Raw and Type1 in 74X samples \n");
400  const std::vector<PackedMETUncertainty> &v = (level == Type1 ? uncertaintiesType1_ : uncertaintiesRaw_);
401  if (v.empty()) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type\n");
402  if (v.size() == 1) {
403  if (shift != MET::METUncertainty::NoShift) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type (only central value available)\n");
404  return Vector2{ (px() + v.front().dpx()), (py() + v.front().dpy()) };
405  }
406  Vector2 ret{ (px() + v[shift].dpx()), (py() + v[shift].dpy()) };
407  return ret;
408 }
409 
411  if (level != Type1 && level != Raw) throw cms::Exception("Unsupported", "MET uncertainties only supported for Raw and Type1 in 74X samples \n");
412  const std::vector<PackedMETUncertainty> &v = (level == Type1 ? uncertaintiesType1_ : uncertaintiesRaw_);
413  if (v.empty()) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type\n");
414  if (v.size() == 1) {
415  if (shift != MET::METUncertainty::NoShift) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type (only central value available)\n");
416  return Vector(px() + v.front().dpx(), py() + v.front().dpy(), 0);
417  }
418  return Vector(px() + v[shift].dpx(), py() + v[shift].dpy(), 0);
419 }
420 
422  if (level != Type1 && level != Raw) throw cms::Exception("Unsupported", "MET uncertainties only supported for Raw and Type1 in 74X samples\n");
423  const std::vector<PackedMETUncertainty> &v = (level == Type1 ? uncertaintiesType1_ : uncertaintiesRaw_);
424  if (v.empty()) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type\n");
425  if (v.size() == 1) {
426  if (shift != MET::METUncertainty::NoShift) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type (only central value available)\n");
427  double x = px() + v.front().dpx(), y = py() + v.front().dpy();
428  return LorentzVector(x, y, 0, std::hypot(x,y));
429  }
430  double x = px() + v[shift].dpx(), y = py() + v[shift].dpy();
431  return LorentzVector(x, y, 0, std::hypot(x,y));
432 }
433 
435  if (level != Type1 && level != Raw) throw cms::Exception("Unsupported", "MET uncertainties only supported for Raw and Type1 in 74X samples\n");
436  const std::vector<PackedMETUncertainty> &v = (level == Type1 ? uncertaintiesType1_ : uncertaintiesRaw_);
437  if (v.empty()) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type\n");
438  if (v.size() == 1) {
439  if (shift != MET::METUncertainty::NoShift) throw cms::Exception("Unsupported", "MET uncertainties not available for the specified correction type (only central value available)\n");
440  return sumEt() + v.front().dsumEt();
441  }
442  return sumEt() + v[shift].dsumEt();
443 }
444 
445 
446 
447 
449 
451  packedDpx_ = MiniFloatConverter::float32to16(dpx_);
452  packedDpy_ = MiniFloatConverter::float32to16(dpy_);
453  packedDSumEt_ = MiniFloatConverter::float32to16(dsumEt_);
454 }
456  unpacked_=true;
457  dpx_=MiniFloatConverter::float16to32(packedDpx_);
458  dpy_=MiniFloatConverter::float16to32(packedDpy_);
459  dsumEt_=MiniFloatConverter::float16to32(packedDSumEt_);
460 }
461 
Analysis-level MET class.
Definition: MET.h:43
std::vector< double > dsumEt() const
Definition: MET.cc:126
value_type const * get() const
Definition: RefToBase.h:234
std::vector< PackedMETUncertainty > uncertaintiesRaw_
Definition: MET.h:265
Vector2 shiftedP2_74x(METUncertainty shift, METCorrectionLevel level) const
Definition: MET.cc:398
const PackedMETUncertainty findMETTotalShift(MET::METCorrectionLevel cor, MET::METUncertainty shift) const
Definition: MET.cc:199
MET()
default constructor
Definition: MET.cc:11
SpecificPFMETData getSpecific() const
Definition: PFMET.h:72
double dpy() const
Definition: MET.h:230
math::XYZVector Vector
point in the space
Definition: Candidate.h:43
double metSig_
Definition: MET.h:252
Vector2 caloMETP2() const
Definition: MET.cc:381
std::vector< PackedMETUncertainty > uncertaintiesType1p2_
Definition: MET.h:265
Vector2 corP2(METCorrectionLevel level=Type1) const
Definition: MET.cc:328
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:159
std::vector< PackedMETUncertainty > uncertainties_
Definition: MET.h:267
double px() const final
x coordinate of momentum vector
void setSignificanceMatrix(const reco::METCovMatrix &matrix)
Definition: MET.cc:157
#define nullptr
double y() const final
rapidity
LorentzVector shiftedP4(METUncertainty shift, METCorrectionLevel level=Type1) const
Definition: MET.cc:281
void setUncShift(double px, double py, double sumEt, METUncertainty shift, bool isSmeared=false)
Definition: MET.cc:355
double dsumEt() const
Definition: MET.h:231
~MET() override
destructor
Definition: MET.cc:89
double metSignificance() const
Definition: MET.cc:126
METUncertainty
Definition: MET.h:151
Vector2 uncorP2() const
Definition: MET.cc:341
Vector corP3(METCorrectionLevel level=Type1) const
Definition: MET.cc:331
static float float16to32(uint16_t h)
Definition: libminifloat.h:12
double caloMETSumEt() const
Definition: MET.cc:393
Definition: HeavyIon.h:7
double corSumEt(METCorrectionLevel level=Type1) const
Definition: MET.cc:337
void add(float dpx, float dpy, float dsumEt)
Definition: MET.h:233
SpecificCaloMETData getSpecific() const
Definition: CaloMET.h:79
double sumEt() const
Definition: MET.h:56
Definition: MET.h:42
PackedMETUncertainty caloPackedMet_
Definition: MET.h:270
static uint16_t float32to16(float x)
Definition: libminifloat.h:17
LorentzVector corP4(METCorrectionLevel level=Type1) const
Definition: MET.cc:334
double dpx() const
Definition: MET.h:229
Vector2 shiftedP2(METUncertainty shift, METCorrectionLevel level=Type1) const
Definition: MET.cc:235
double caloMETPhi() const
Definition: MET.cc:389
double shiftedSumEt(METUncertainty shift, METCorrectionLevel level=Type1) const
Definition: MET.cc:306
double phi() const
Definition: MET.h:172
std::vector< reco::GenMET > genMET_
Definition: MET.h:245
LorentzVector shiftedP4_74x(METUncertainty shift, METCorrectionLevel level) const
Definition: MET.cc:421
void initCorMap()
Definition: MET.cc:131
double shiftedSumEt_74x(METUncertainty shift, METCorrectionLevel level) const
Definition: MET.cc:434
double uncorSumEt() const
Definition: MET.cc:350
double pt() const
Definition: MET.h:171
MET & operator=(MET const &)
Definition: MET.cc:93
std::map< MET::METCorrectionLevel, std::vector< MET::METCorrectionType > > corMap_
Definition: MET.h:256
const reco::GenMET * genMET() const
return the associated GenMET
Definition: MET.cc:110
double py() const final
y coordinate of momentum vector
std::vector< SpecificCaloMETData > caloMET_
Definition: MET.h:247
this below should be private but Reflex doesn&#39;t like it
Definition: MET.h:223
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
METCorrectionLevel
Definition: MET.h:158
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
std::vector< PackedMETUncertainty > uncertaintiesType1_
Definition: MET.h:265
fixed size matrix
METCorrectionType
Definition: MET.h:163
std::vector< PackedMETUncertainty > corrections_
Definition: MET.h:268
math::XYZVector Vector
point in the space
Definition: LeafCandidate.h:29
static unsigned int const shift
LorentzVector uncorP4() const
Definition: MET.cc:347
Templated PAT object container.
Definition: PATObject.h:49
void setCorShift(double px, double py, double sumEt, METCorrectionType level)
Definition: MET.cc:372
void setGenMET(const reco::GenMET &gm)
set the associated GenMET
Definition: MET.cc:115
std::vector< SpecificPFMETData > pfMET_
Definition: MET.h:249
double caloMETPt() const
Definition: MET.cc:385
Vector uncorP3() const
Definition: MET.cc:344
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
void unpack() const
Definition: MET.cc:455
void setMETSignificance(const double &metSig)
Definition: MET.cc:122
Vector shiftedP3_74x(METUncertainty shift, METCorrectionLevel level) const
Definition: MET.cc:410
reco::METCovMatrix getSignificanceMatrix(void) const
Definition: MET.cc:139
Vector shiftedP3(METUncertainty shift, METCorrectionLevel level=Type1) const
Definition: MET.cc:258