00001 #include "DataFormats/METReco/interface/MET.h" 00002 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00003 #include "DQMOffline/PFTau/interface/Matchers.h" 00004 00005 #include "DQMOffline/PFTau/interface/PFMETMonitor.h" 00006 00007 #include <TROOT.h> 00008 #include <TFile.h> 00009 #include <TH1.h> 00010 #include <TH2.h> 00011 00012 PFMETMonitor::~PFMETMonitor() {} 00013 00014 00015 void PFMETMonitor::setParameters( const edm::ParameterSet & parameterSet) { 00016 00017 mode_ = (Benchmark::Mode) parameterSet.getParameter<int>( "mode" ); 00018 variablePtBins_ = parameterSet.getParameter< std::vector<double> >( "VariablePtBins" ); 00019 createMETSpecificHistos_ = parameterSet.getParameter<bool>( "CreateMETSpecificHistos" ); 00020 setRange( parameterSet.getParameter<double>("ptMin"), 00021 parameterSet.getParameter<double>("ptMax"), 00022 parameterSet.getParameter<double>("etaMin"), 00023 parameterSet.getParameter<double>("etaMax"), 00024 parameterSet.getParameter<double>("phiMin"), 00025 parameterSet.getParameter<double>("phiMax") ); 00026 00027 00028 candBench_.setParameters(mode_); 00029 matchCandBench_.setParameters(mode_); 00030 00031 px_ = 0; 00032 sumEt_ = 0; 00033 delta_ex_ = 0; 00034 delta_ex_VS_set_ = 0; 00035 delta_set_VS_set_ = 0; 00036 delta_set_Over_set_VS_set_ = 0; 00037 } 00038 00039 void PFMETMonitor::setup() { 00040 candBench_.setup(); 00041 matchCandBench_.setup(); 00042 00043 if (createMETSpecificHistos_) { 00044 float* ptBins = new float[variablePtBins_.size()]; 00045 for (size_t i = 0; i < variablePtBins_.size(); i++) { 00046 ptBins[i] = variablePtBins_[i]; 00047 } 00048 PhaseSpace pxPS = PhaseSpace( 50, 0, 200); 00049 PhaseSpace dpxPS = PhaseSpace( 50, -500, 500); 00050 PhaseSpace dptPS = PhaseSpace( 200, -500, 500); 00051 PhaseSpace setPS = PhaseSpace( 50, 0.0, 3000); 00052 PhaseSpace dsetPS = PhaseSpace( 50, -1000.0, 1000); 00053 PhaseSpace setOvsetPS = PhaseSpace( 100,0., 2.); 00054 00055 px_ = book1D("px_", "px_;p_{X} (GeV)", pxPS.n, pxPS.m, pxPS.M); 00056 00057 sumEt_ = book1D("sumEt_", "sumEt_;#sumE_{T}", setPS.n, setPS.m, setPS.M); 00058 00059 delta_ex_ = book1D("delta_ex_", "#DeltaME_{X}", dpxPS.n, dpxPS.m, dpxPS.M); 00060 00061 delta_ex_VS_set_ = book2D("delta_ex_VS_set_", ";SE_{T, true} (GeV);#DeltaE_{X}", 00062 setPS.n, setPS.m, setPS.M, 00063 dpxPS.n, dpxPS.m, dpxPS.M ); 00064 00065 delta_set_VS_set_ = book2D("delta_set_VS_set_", 00066 ";SE_{T, true} (GeV);#DeltaSE_{T}", 00067 setPS.n, setPS.m, setPS.M, 00068 dsetPS.n, dsetPS.m, dsetPS.M ); 00069 00070 delta_set_Over_set_VS_set_ = book2D("delta_set_Over_set_VS_set_", 00071 ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}", 00072 setPS.n, setPS.m, setPS.M, 00073 setOvsetPS.n, setOvsetPS.m, setOvsetPS.M ); 00074 } 00075 } 00076 void PFMETMonitor::setup(const edm::ParameterSet & parameterSet) { 00077 candBench_.setup(parameterSet); 00078 matchCandBench_.setup(parameterSet); 00079 00080 if (createMETSpecificHistos_) { 00081 float* ptBins = new float[variablePtBins_.size()]; 00082 for (size_t i = 0; i < variablePtBins_.size(); i++) { 00083 ptBins[i] = variablePtBins_[i]; 00084 } 00085 00086 edm::ParameterSet pxPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPxHistoParameter"); 00087 edm::ParameterSet dpxPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPxHistoParameter"); 00088 edm::ParameterSet dptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtHistoParameter"); 00089 edm::ParameterSet setPS = parameterSet.getParameter<edm::ParameterSet>("SumEtHistoParameter"); 00090 edm::ParameterSet dsetPS = parameterSet.getParameter<edm::ParameterSet>("DeltaSumEtHistoParameter"); 00091 edm::ParameterSet setOvsetPS = parameterSet.getParameter<edm::ParameterSet>("DeltaSumEtOvSumEtHistoParameter"); 00092 00093 if (pxPS.getParameter<bool>("switchOn")) { 00094 px_ = book1D("px_", "px_;p_{X} (GeV)", 00095 pxPS.getParameter<int32_t>("nBin"), 00096 pxPS.getParameter<double>("xMin"), 00097 pxPS.getParameter<double>("xMax")); 00098 } 00099 if (setPS.getParameter<bool>("switchOn")) { 00100 sumEt_ = book1D("sumEt_", "sumEt_;#sumE_{T}", 00101 setPS.getParameter<int32_t>("nBin"), 00102 setPS.getParameter<double>("xMin"), 00103 setPS.getParameter<double>("xMax")); 00104 } 00105 if (dpxPS.getParameter<bool>("switchOn")) { 00106 delta_ex_ = book1D("delta_ex_", "#DeltaME_{X}", 00107 dpxPS.getParameter<int32_t>("nBin"), 00108 dpxPS.getParameter<double>("xMin"), 00109 dpxPS.getParameter<double>("xMax")); 00110 } 00111 00112 if (dpxPS.getParameter<bool>("switchOn")) { 00113 delta_ex_VS_set_ = book2D("delta_ex_VS_set_", 00114 ";SE_{T, true} (GeV);#DeltaE_{X}", 00115 setPS.getParameter<int32_t>("nBin"), 00116 setPS.getParameter<double>("xMin"), 00117 setPS.getParameter<double>("xMax"), 00118 dptPS.getParameter<int32_t>("nBin"), 00119 dptPS.getParameter<double>("xMin"), 00120 dptPS.getParameter<double>("xMax")); 00121 } 00122 if (dsetPS.getParameter<bool>("switchOn")) { 00123 delta_set_VS_set_ = book2D("delta_set_VS_set_", 00124 ";SE_{T, true} (GeV);#DeltaSE_{T}", 00125 setPS.getParameter<int32_t>("nBin"), 00126 setPS.getParameter<double>("xMin"), 00127 setPS.getParameter<double>("xMax"), 00128 dsetPS.getParameter<int32_t>("nBin"), 00129 dsetPS.getParameter<double>("xMin"), 00130 dsetPS.getParameter<double>("xMax")); 00131 } 00132 if (setOvsetPS.getParameter<bool>("switchOn")) { 00133 delta_set_Over_set_VS_set_ = book2D("delta_set_Over_set_VS_set_", 00134 ";SE_{T, true} (GeV);#DeltaSE_{T}/SE_{T}", 00135 setPS.getParameter<int32_t>("nBin"), 00136 setPS.getParameter<double>("xMin"), 00137 setPS.getParameter<double>("xMax"), 00138 setOvsetPS.getParameter<int32_t>("nBin"), 00139 setOvsetPS.getParameter<double>("xMin"), 00140 setOvsetPS.getParameter<double>("xMax")); 00141 } 00142 } 00143 } 00144 00145 void PFMETMonitor::setDirectory(TDirectory* dir) { 00146 Benchmark::setDirectory(dir); 00147 00148 candBench_.setDirectory(dir); 00149 matchCandBench_.setDirectory(dir); 00150 } 00151 00152 void PFMETMonitor::fillOne(const reco::MET& met, 00153 const reco::MET& matchedMet, float& minVal, float& maxVal) { 00154 candBench_.fillOne(met); 00155 matchCandBench_.fillOne(met, matchedMet); 00156 if (createMETSpecificHistos_) { 00157 if( !isInRange(met.pt(), met.eta(), met.phi() ) ) return; 00158 00159 if (px_) px_->Fill(met.px()); 00160 if (delta_ex_) { 00161 delta_ex_->Fill(met.px()-matchedMet.px()); 00162 delta_ex_->Fill(met.py()-matchedMet.py()); 00163 } 00164 if (sumEt_) sumEt_->Fill( met.sumEt()); 00165 00166 if (delta_ex_VS_set_) { 00167 delta_ex_VS_set_->Fill(matchedMet.sumEt(), met.px()-matchedMet.px()); 00168 delta_ex_VS_set_->Fill(matchedMet.sumEt(), met.py()-matchedMet.py()); 00169 } 00170 if (delta_set_VS_set_) delta_set_VS_set_->Fill(matchedMet.sumEt(),met.sumEt()-matchedMet.sumEt()); 00171 if (delta_set_Over_set_VS_set_ && matchedMet.sumEt()>0.001 ) { 00172 float setRes = (met.sumEt() - matchedMet.sumEt())/matchedMet.sumEt(); 00173 if (setRes > maxVal) maxVal = setRes; 00174 if (setRes < minVal) minVal = setRes; 00175 delta_set_Over_set_VS_set_->Fill(matchedMet.sumEt(),setRes); 00176 } 00177 } 00178 }