Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:58:19

0001 #include "DQMOffline/PFTau/interface/MatchMETBenchmark.h"
0002 #include "DataFormats/Candidate/interface/Candidate.h"
0003 #include "DataFormats/METReco/interface/MET.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 
0006 #include <TFile.h>
0007 #include <TH1.h>
0008 #include <TH2.h>
0009 #include <TROOT.h>
0010 
0011 using namespace std;
0012 
0013 MatchMETBenchmark::~MatchMETBenchmark() {}
0014 
0015 void MatchMETBenchmark::setup(DQMStore::IBooker &b) {
0016   // std::cout << "FL: MatchMETBenchmark.cc: start setup()" << std::endl;
0017   PhaseSpace ptPS;
0018   PhaseSpace dptOvptPS;
0019   PhaseSpace dptPS;
0020   PhaseSpace dphiPS;
0021   PhaseSpace setPS;
0022   PhaseSpace dsetPS;
0023   PhaseSpace setOvsetPS;
0024 
0025   switch (mode_) {
0026     case VALIDATION:
0027       ptPS = PhaseSpace(100, 0, 1000);
0028       dptOvptPS = PhaseSpace(200, -1, 1);
0029       dphiPS = PhaseSpace(100, -3.2, 3.2);
0030       dptPS = PhaseSpace(200, -500, 500);
0031       setPS = PhaseSpace(300, 0.0, 3000);
0032       dsetPS = PhaseSpace(200, 0. - 1000, 1000);
0033       setOvsetPS = PhaseSpace(500, 0., 2.);
0034       break;
0035     case DQMOFFLINE:
0036     default:
0037       ptPS = PhaseSpace(50, 0, 200);
0038       dptOvptPS = PhaseSpace(50, -1, 1);
0039       dphiPS = PhaseSpace(50, -3.2, 3.2);
0040       dptPS = PhaseSpace(50, -500, 500);
0041       setPS = PhaseSpace(50, 0.0, 3000);
0042       dsetPS = PhaseSpace(50, -1000.0, 1000);
0043       setOvsetPS = PhaseSpace(100, 0., 2.);
0044       break;
0045   }
0046 
0047   // variable bins to be done here, as they will save a lot of memory.
0048 
0049   // float ptBins[11] = {0, 1, 2, 5, 10, 20, 50, 100, 200, 400, 1000};
0050 
0051   delta_et_Over_et_VS_et_ = book2D(b,
0052                                    "delta_et_Over_et_VS_et_",
0053                                    ";ME_{T, true} (GeV);#DeltaME_{T}/ME_{T}",
0054                                    ptPS.n,
0055                                    ptPS.m,
0056                                    ptPS.M,
0057                                    dptOvptPS.n,
0058                                    dptOvptPS.m,
0059                                    dptOvptPS.M);
0060 
0061   delta_et_VS_et_ = book2D(
0062       b, "delta_et_VS_et_", ";ME_{T, true} (GeV);#DeltaME_{T}", ptPS.n, ptPS.m, ptPS.M, dptPS.n, dptPS.m, dptPS.M);
0063 
0064   delta_phi_VS_et_ = book2D(
0065       b, "delta_phi_VS_et_", ";ME_{T, true} (GeV);#Delta#phi", ptPS.n, ptPS.m, ptPS.M, dphiPS.n, dphiPS.m, dphiPS.M);
0066 
0067   delta_ex_ = book1D(b, "delta_ex_", "#DeltaME_{X}", dptPS.n, dptPS.m, dptPS.M);
0068 
0069   RecEt_VS_TrueEt_ =
0070       book2D(b, "RecEt_VS_TrueEt_", ";ME_{T, true} (GeV);ME_{T}", ptPS.n, ptPS.m, ptPS.M, ptPS.n, ptPS.m, ptPS.M);
0071 
0072   delta_set_VS_set_ = book2D(b,
0073                              "delta_set_VS_set_",
0074                              ";SE_{T, true} (GeV);#DeltaSE_{T}",
0075                              setPS.n,
0076                              setPS.m,
0077                              setPS.M,
0078                              dsetPS.n,
0079                              dsetPS.m,
0080                              dsetPS.M);
0081 
0082   delta_set_Over_set_VS_set_ = book2D(b,
0083                                       "delta_set_Over_set_VS_set_",
0084                                       ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}",
0085                                       setPS.n,
0086                                       setPS.m,
0087                                       setPS.M,
0088                                       dptOvptPS.n,
0089                                       dptOvptPS.m,
0090                                       dptOvptPS.M);
0091 
0092   delta_ex_VS_set_ = book2D(
0093       b, "delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}", setPS.n, setPS.m, setPS.M, ptPS.n, -ptPS.M, ptPS.M);
0094 
0095   RecSet_Over_TrueSet_VS_TrueSet_ = book2D(b,
0096                                            "RecSet_Over_TrueSet_VS_TrueSet_",
0097                                            ";SE_{T, true} (GeV);SE_{T}/SE_{T}",
0098                                            setPS.n,
0099                                            setPS.m,
0100                                            setPS.M,
0101                                            setOvsetPS.n,
0102                                            setOvsetPS.m,
0103                                            setOvsetPS.M);
0104 }
0105 
0106 void MatchMETBenchmark::fillOne(const reco::MET &cand, const reco::MET &matchedCand) {
0107   if (!isInRange(cand.pt(), cand.eta(), cand.phi()))
0108     return;
0109 
0110   if (matchedCand.pt() > 0.001)
0111     delta_et_Over_et_VS_et_->Fill(matchedCand.pt(), (cand.pt() - matchedCand.pt()) / matchedCand.pt());
0112   else
0113     edm::LogWarning("MatchMETBenchmark") << " matchedCand.pt()<0.001";
0114   delta_et_VS_et_->Fill(matchedCand.pt(), cand.pt() - matchedCand.pt());
0115   delta_phi_VS_et_->Fill(matchedCand.pt(), cand.phi() - matchedCand.phi());
0116   delta_ex_->Fill(cand.px() - matchedCand.px());
0117   delta_ex_->Fill(cand.py() - matchedCand.py());
0118   RecEt_VS_TrueEt_->Fill(matchedCand.pt(), cand.pt());
0119   delta_set_VS_set_->Fill(matchedCand.sumEt(), cand.sumEt() - matchedCand.sumEt());
0120   if (matchedCand.sumEt() > 0.001)
0121     delta_set_Over_set_VS_set_->Fill(matchedCand.sumEt(), (cand.sumEt() - matchedCand.sumEt()) / matchedCand.sumEt());
0122   else
0123     edm::LogWarning("MatchMETBenchmark") << " matchedCand.sumEt()<0.001";
0124   delta_ex_VS_set_->Fill(matchedCand.sumEt(), cand.px() - matchedCand.px());
0125   delta_ex_VS_set_->Fill(matchedCand.sumEt(), cand.py() - matchedCand.py());
0126   if (matchedCand.sumEt() > 0.001)
0127     RecSet_Over_TrueSet_VS_TrueSet_->Fill(matchedCand.sumEt(), cand.sumEt() / matchedCand.sumEt());
0128 }