File indexing completed on 2024-04-06 12:33:39
0001 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0002 #include "DQMServices/Core/interface/DQMStore.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009
0010 class Phase2OTHarvestTrackingParticles : public DQMEDHarvester {
0011 public:
0012 explicit Phase2OTHarvestTrackingParticles(const edm::ParameterSet &);
0013 ~Phase2OTHarvestTrackingParticles() override;
0014 void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override;
0015 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0016
0017 private:
0018
0019 DQMStore *dbe;
0020 std::string topFolderName_;
0021 };
0022
0023 Phase2OTHarvestTrackingParticles::Phase2OTHarvestTrackingParticles(const edm::ParameterSet &iConfig)
0024 : topFolderName_(iConfig.getParameter<std::string>("TopFolderName")) {}
0025
0026 Phase2OTHarvestTrackingParticles::~Phase2OTHarvestTrackingParticles() {}
0027
0028
0029
0030 void Phase2OTHarvestTrackingParticles::dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0031 using namespace edm;
0032
0033 float eta_bins[] = {0.0, 0.7, 1.0, 1.2, 1.6, 2.0, 2.4};
0034 int eta_binnum = 6;
0035
0036 dbe = nullptr;
0037 dbe = edm::Service<DQMStore>().operator->();
0038
0039 if (dbe) {
0040
0041 MonitorElement *meN_eta = dbe->get(topFolderName_ + "/Efficiency/match_tp_eta");
0042 MonitorElement *meD_eta = dbe->get(topFolderName_ + "/Efficiency/tp_eta");
0043 MonitorElement *meN_pt = dbe->get(topFolderName_ + "/Efficiency/match_tp_pt");
0044 MonitorElement *meD_pt = dbe->get(topFolderName_ + "/Efficiency/tp_pt");
0045 MonitorElement *meN_pt_zoom = dbe->get(topFolderName_ + "/Efficiency/match_tp_pt_zoom");
0046 MonitorElement *meD_pt_zoom = dbe->get(topFolderName_ + "/Efficiency/tp_pt_zoom");
0047 MonitorElement *meN_d0 = dbe->get(topFolderName_ + "/Efficiency/match_tp_d0");
0048 MonitorElement *meD_d0 = dbe->get(topFolderName_ + "/Efficiency/tp_d0");
0049 MonitorElement *meN_VtxR = dbe->get(topFolderName_ + "/Efficiency/match_tp_VtxR");
0050 MonitorElement *meD_VtxR = dbe->get(topFolderName_ + "/Efficiency/tp_VtxR");
0051 MonitorElement *meN_VtxZ = dbe->get(topFolderName_ + "/Efficiency/match_tp_VtxZ");
0052 MonitorElement *meD_VtxZ = dbe->get(topFolderName_ + "/Efficiency/tp_VtxZ");
0053
0054 MonitorElement *merespt_eta0to0p7_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta0to0p7_pt2to3");
0055 MonitorElement *merespt_eta0p7to1_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta0p7to1_pt2to3");
0056 MonitorElement *merespt_eta1to1p2_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta1to1p2_pt2to3");
0057 MonitorElement *merespt_eta1p2to1p6_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p2to1p6_pt2to3");
0058 MonitorElement *merespt_eta1p6to2_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p6to2_pt2to3");
0059 MonitorElement *merespt_eta2to2p4_pt2to3 = dbe->get(topFolderName_ + "/Resolution/respt_eta2to2p4_pt2to3");
0060 MonitorElement *merespt_eta0to0p7_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta0to0p7_pt3to8");
0061 MonitorElement *merespt_eta0p7to1_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta0p7to1_pt3to8");
0062 MonitorElement *merespt_eta1to1p2_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta1to1p2_pt3to8");
0063 MonitorElement *merespt_eta1p2to1p6_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p2to1p6_pt3to8");
0064 MonitorElement *merespt_eta1p6to2_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta1p6to2_pt3to8");
0065 MonitorElement *merespt_eta2to2p4_pt3to8 = dbe->get(topFolderName_ + "/Resolution/respt_eta2to2p4_pt3to8");
0066 MonitorElement *merespt_eta0to0p7_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta0to0p7_pt8toInf");
0067 MonitorElement *merespt_eta0p7to1_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta0p7to1_pt8toInf");
0068 MonitorElement *merespt_eta1to1p2_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta1to1p2_pt8toInf");
0069 MonitorElement *merespt_eta1p2to1p6_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta1p2to1p6_pt8toInf");
0070 MonitorElement *merespt_eta1p6to2_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta1p6to2_pt8toInf");
0071 MonitorElement *merespt_eta2to2p4_pt8toInf = dbe->get(topFolderName_ + "/Resolution/respt_eta2to2p4_pt8toInf");
0072
0073 MonitorElement *mereseta_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/reseta_eta0to0p7");
0074 MonitorElement *mereseta_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/reseta_eta0p7to1");
0075 MonitorElement *mereseta_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/reseta_eta1to1p2");
0076 MonitorElement *mereseta_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/reseta_eta1p2to1p6");
0077 MonitorElement *mereseta_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/reseta_eta1p6to2");
0078 MonitorElement *mereseta_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/reseta_eta2to2p4");
0079
0080 MonitorElement *meresphi_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/resphi_eta0to0p7");
0081 MonitorElement *meresphi_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/resphi_eta0p7to1");
0082 MonitorElement *meresphi_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/resphi_eta1to1p2");
0083 MonitorElement *meresphi_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/resphi_eta1p2to1p6");
0084 MonitorElement *meresphi_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/resphi_eta1p6to2");
0085 MonitorElement *meresphi_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/resphi_eta2to2p4");
0086
0087 MonitorElement *meresVtxZ_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta0to0p7");
0088 MonitorElement *meresVtxZ_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta0p7to1");
0089 MonitorElement *meresVtxZ_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta1to1p2");
0090 MonitorElement *meresVtxZ_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta1p2to1p6");
0091 MonitorElement *meresVtxZ_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta1p6to2");
0092 MonitorElement *meresVtxZ_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/resVtxZ_eta2to2p4");
0093
0094 MonitorElement *meresd0_eta0to0p7 = dbe->get(topFolderName_ + "/Resolution/resd0_eta0to0p7");
0095 MonitorElement *meresd0_eta0p7to1 = dbe->get(topFolderName_ + "/Resolution/resd0_eta0p7to1");
0096 MonitorElement *meresd0_eta1to1p2 = dbe->get(topFolderName_ + "/Resolution/resd0_eta1to1p2");
0097 MonitorElement *meresd0_eta1p2to1p6 = dbe->get(topFolderName_ + "/Resolution/resd0_eta1p2to1p6");
0098 MonitorElement *meresd0_eta1p6to2 = dbe->get(topFolderName_ + "/Resolution/resd0_eta1p6to2");
0099 MonitorElement *meresd0_eta2to2p4 = dbe->get(topFolderName_ + "/Resolution/resd0_eta2to2p4");
0100
0101 if (meN_eta && meD_eta) {
0102
0103 TH1F *numerator = meN_eta->getTH1F();
0104 TH1F *denominator = meD_eta->getTH1F();
0105 numerator->Sumw2();
0106 denominator->Sumw2();
0107
0108
0109 dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0110
0111
0112 MonitorElement *me_effic_eta = ibooker.book1D("EtaEfficiency",
0113 "#eta efficiency",
0114 numerator->GetNbinsX(),
0115 numerator->GetXaxis()->GetXmin(),
0116 numerator->GetXaxis()->GetXmax());
0117
0118
0119 me_effic_eta->getTH1F()->Divide(numerator, denominator, 1., 1., "B");
0120 me_effic_eta->setAxisTitle("tracking particle #eta");
0121 me_effic_eta->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0122 me_effic_eta->getTH1F()->SetMaximum(1.0);
0123 me_effic_eta->getTH1F()->SetMinimum(0.0);
0124 me_effic_eta->getTH1F()->SetStats(false);
0125 }
0126 else {
0127 edm::LogWarning("DataNotFound") << "Monitor elements for eta efficiency cannot be found!\n";
0128 }
0129
0130 if (meN_pt && meD_pt) {
0131
0132 TH1F *numerator2 = meN_pt->getTH1F();
0133 numerator2->Sumw2();
0134 TH1F *denominator2 = meD_pt->getTH1F();
0135 denominator2->Sumw2();
0136
0137
0138 dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0139
0140
0141 MonitorElement *me_effic_pt = ibooker.book1D("PtEfficiency",
0142 "p_{T} efficiency",
0143 numerator2->GetNbinsX(),
0144 numerator2->GetXaxis()->GetXmin(),
0145 numerator2->GetXaxis()->GetXmax());
0146
0147
0148 me_effic_pt->getTH1F()->Divide(numerator2, denominator2, 1., 1., "B");
0149 me_effic_pt->setAxisTitle("Tracking particle p_{T} [GeV]");
0150 me_effic_pt->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0151 me_effic_pt->getTH1F()->SetMaximum(1.0);
0152 me_effic_pt->getTH1F()->SetMinimum(0.0);
0153 me_effic_pt->getTH1F()->SetStats(false);
0154 }
0155 else {
0156 edm::LogWarning("DataNotFound") << "Monitor elements for pT efficiency cannot be found!\n";
0157 }
0158
0159 if (meN_pt_zoom && meD_pt_zoom) {
0160
0161 TH1F *numerator2_zoom = meN_pt_zoom->getTH1F();
0162 numerator2_zoom->Sumw2();
0163 TH1F *denominator2_zoom = meD_pt_zoom->getTH1F();
0164 denominator2_zoom->Sumw2();
0165
0166
0167 dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0168
0169
0170 MonitorElement *me_effic_pt_zoom = ibooker.book1D("PtEfficiency_zoom",
0171 "p_{T} efficiency",
0172 numerator2_zoom->GetNbinsX(),
0173 numerator2_zoom->GetXaxis()->GetXmin(),
0174 numerator2_zoom->GetXaxis()->GetXmax());
0175
0176
0177 me_effic_pt_zoom->getTH1F()->Divide(numerator2_zoom, denominator2_zoom, 1., 1., "B");
0178 me_effic_pt_zoom->setAxisTitle("Tracking particle p_{T} [GeV]");
0179 me_effic_pt_zoom->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0180 me_effic_pt_zoom->getTH1F()->SetMaximum(1.0);
0181 me_effic_pt_zoom->getTH1F()->SetMinimum(0.0);
0182 me_effic_pt_zoom->getTH1F()->SetStats(false);
0183 }
0184 else {
0185 edm::LogWarning("DataNotFound") << "Monitor elements for zoom pT efficiency cannot be found!\n";
0186 }
0187
0188 if (meN_d0 && meD_d0) {
0189
0190 TH1F *numerator5 = meN_d0->getTH1F();
0191 numerator5->Sumw2();
0192 TH1F *denominator5 = meD_d0->getTH1F();
0193 denominator5->Sumw2();
0194
0195
0196 dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0197
0198
0199 MonitorElement *me_effic_d0 = ibooker.book1D("d0Efficiency",
0200 "d_{0} efficiency",
0201 numerator5->GetNbinsX(),
0202 numerator5->GetXaxis()->GetXmin(),
0203 numerator5->GetXaxis()->GetXmax());
0204
0205
0206 me_effic_d0->getTH1F()->Divide(numerator5, denominator5, 1., 1., "B");
0207 me_effic_d0->setAxisTitle("Tracking particle d_{0} [cm]");
0208 me_effic_d0->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0209 me_effic_d0->getTH1F()->SetMaximum(1.0);
0210 me_effic_d0->getTH1F()->SetMinimum(0.0);
0211 me_effic_d0->getTH1F()->SetStats(false);
0212 }
0213 else {
0214 edm::LogWarning("DataNotFound") << "Monitor elements for d0 efficiency cannot be found!\n";
0215 }
0216
0217 if (meN_VtxR && meD_VtxR) {
0218
0219 TH1F *numerator6 = meN_VtxR->getTH1F();
0220 numerator6->Sumw2();
0221 TH1F *denominator6 = meD_VtxR->getTH1F();
0222 denominator6->Sumw2();
0223
0224
0225 dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0226
0227
0228 MonitorElement *me_effic_VtxR = ibooker.book1D("VtxREfficiency",
0229 "Vtx R efficiency",
0230 numerator6->GetNbinsX(),
0231 numerator6->GetXaxis()->GetXmin(),
0232 numerator6->GetXaxis()->GetXmax());
0233
0234
0235 me_effic_VtxR->getTH1F()->Divide(numerator6, denominator6, 1., 1., "B");
0236 me_effic_VtxR->setAxisTitle("Tracking particle VtxR [cm]");
0237 me_effic_VtxR->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0238 me_effic_VtxR->getTH1F()->SetMaximum(1.0);
0239 me_effic_VtxR->getTH1F()->SetMinimum(0.0);
0240 me_effic_VtxR->getTH1F()->SetStats(false);
0241 }
0242 else {
0243 edm::LogWarning("DataNotFound") << "Monitor elements for VtxR efficiency cannot be found!\n";
0244 }
0245
0246 if (meN_VtxZ && meD_VtxZ) {
0247
0248 TH1F *numerator7 = meN_VtxZ->getTH1F();
0249 numerator7->Sumw2();
0250 TH1F *denominator7 = meD_VtxZ->getTH1F();
0251 denominator7->Sumw2();
0252
0253
0254 dbe->setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0255
0256
0257 MonitorElement *me_effic_VtxZ = ibooker.book1D("VtxZEfficiency",
0258 "Vtx Z efficiency",
0259 numerator7->GetNbinsX(),
0260 numerator7->GetXaxis()->GetXmin(),
0261 numerator7->GetXaxis()->GetXmax());
0262
0263
0264 me_effic_VtxZ->getTH1F()->Divide(numerator7, denominator7, 1., 1., "B");
0265 me_effic_VtxZ->setAxisTitle("Tracking particle VtxZ [cm]");
0266 me_effic_VtxZ->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0267 me_effic_VtxZ->getTH1F()->SetMaximum(1.0);
0268 me_effic_VtxZ->getTH1F()->SetMinimum(0.0);
0269 me_effic_VtxZ->getTH1F()->SetStats(false);
0270 }
0271 else {
0272 edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ efficiency cannot be found!\n";
0273 }
0274
0275 if (merespt_eta0to0p7_pt2to3 && merespt_eta0p7to1_pt2to3 && merespt_eta1to1p2_pt2to3 &&
0276 merespt_eta1p2to1p6_pt2to3 && merespt_eta1p6to2_pt2to3 && merespt_eta2to2p4_pt2to3) {
0277
0278 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0279
0280
0281 TH1F *resPt1a = merespt_eta0to0p7_pt2to3->getTH1F();
0282 TH1F *resPt2a = merespt_eta0p7to1_pt2to3->getTH1F();
0283 TH1F *resPt3a = merespt_eta1to1p2_pt2to3->getTH1F();
0284 TH1F *resPt4a = merespt_eta1p2to1p6_pt2to3->getTH1F();
0285 TH1F *resPt5a = merespt_eta1p6to2_pt2to3->getTH1F();
0286 TH1F *resPt6a = merespt_eta2to2p4_pt2to3->getTH1F();
0287
0288
0289 MonitorElement *me_res_pt1 =
0290 ibooker.book1D("pTResVsEta_2-3", "p_{T} resolution vs |#eta|, for p_{T}: 2-3 GeV", eta_binnum, eta_bins);
0291 TH1F *resPt1 = me_res_pt1->getTH1F();
0292 resPt1->GetXaxis()->SetTitle("tracking particle |#eta|");
0293 resPt1->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
0294 resPt1->SetMinimum(0.0);
0295 resPt1->SetStats(false);
0296
0297 std::vector<TH1F *> vResPt1 = {resPt1a, resPt2a, resPt3a, resPt4a, resPt5a, resPt6a};
0298 for (int i = 0; i < 6; i++) {
0299 resPt1->SetBinContent(i + 1, vResPt1[i]->GetStdDev());
0300 resPt1->SetBinError(i + 1, vResPt1[i]->GetStdDevError());
0301 }
0302 }
0303 else {
0304 edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (2-3) cannot be found!\n";
0305 }
0306
0307 if (merespt_eta0to0p7_pt3to8 && merespt_eta0p7to1_pt3to8 && merespt_eta1to1p2_pt3to8 &&
0308 merespt_eta1p2to1p6_pt3to8 && merespt_eta1p6to2_pt3to8 && merespt_eta2to2p4_pt3to8) {
0309
0310 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0311
0312
0313 TH1F *resPt1b = merespt_eta0to0p7_pt3to8->getTH1F();
0314 TH1F *resPt2b = merespt_eta0p7to1_pt3to8->getTH1F();
0315 TH1F *resPt3b = merespt_eta1to1p2_pt3to8->getTH1F();
0316 TH1F *resPt4b = merespt_eta1p2to1p6_pt3to8->getTH1F();
0317 TH1F *resPt5b = merespt_eta1p6to2_pt3to8->getTH1F();
0318 TH1F *resPt6b = merespt_eta2to2p4_pt3to8->getTH1F();
0319
0320
0321 MonitorElement *me_res_pt2 =
0322 ibooker.book1D("pTResVsEta_3-8", "p_{T} resolution vs |#eta|, for p_{T}: 3-8 GeV", eta_binnum, eta_bins);
0323 TH1F *resPt2 = me_res_pt2->getTH1F();
0324 resPt2->GetXaxis()->SetTitle("tracking particle |#eta|");
0325 resPt2->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
0326 resPt2->SetMinimum(0.0);
0327 resPt2->SetStats(false);
0328
0329 std::vector<TH1F *> vResPt2 = {resPt1b, resPt2b, resPt3b, resPt4b, resPt5b, resPt6b};
0330 for (int i = 0; i < 6; i++) {
0331 resPt2->SetBinContent(i + 1, vResPt2[i]->GetStdDev());
0332 resPt2->SetBinError(i + 1, vResPt2[i]->GetStdDevError());
0333 }
0334 }
0335 else {
0336 edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (3-8) cannot be found!\n";
0337 }
0338
0339 if (merespt_eta0to0p7_pt8toInf && merespt_eta0p7to1_pt8toInf && merespt_eta1to1p2_pt8toInf &&
0340 merespt_eta1p2to1p6_pt8toInf && merespt_eta1p6to2_pt8toInf && merespt_eta2to2p4_pt8toInf) {
0341
0342 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0343
0344
0345 TH1F *resPt1c = merespt_eta0to0p7_pt8toInf->getTH1F();
0346 TH1F *resPt2c = merespt_eta0p7to1_pt8toInf->getTH1F();
0347 TH1F *resPt3c = merespt_eta1to1p2_pt8toInf->getTH1F();
0348 TH1F *resPt4c = merespt_eta1p2to1p6_pt8toInf->getTH1F();
0349 TH1F *resPt5c = merespt_eta1p6to2_pt8toInf->getTH1F();
0350 TH1F *resPt6c = merespt_eta2to2p4_pt8toInf->getTH1F();
0351
0352
0353 MonitorElement *me_res_pt3 =
0354 ibooker.book1D("pTResVsEta_8-inf", "p_{T} resolution vs |#eta|, for p_{T}: >8 GeV", eta_binnum, eta_bins);
0355 TH1F *resPt3 = me_res_pt3->getTH1F();
0356 resPt3->GetXaxis()->SetTitle("tracking particle |#eta|");
0357 resPt3->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
0358 resPt3->SetMinimum(0.0);
0359 resPt3->SetStats(false);
0360
0361 std::vector<TH1F *> vResPt3 = {resPt1c, resPt2c, resPt3c, resPt4c, resPt5c, resPt6c};
0362 for (int i = 0; i < 6; i++) {
0363 resPt3->SetBinContent(i + 1, vResPt3[i]->GetStdDev());
0364 resPt3->SetBinError(i + 1, vResPt3[i]->GetStdDevError());
0365 }
0366 }
0367 else {
0368 edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (8-inf) cannot be found!\n";
0369 }
0370
0371 if (mereseta_eta0to0p7 && mereseta_eta0p7to1 && mereseta_eta1to1p2 && mereseta_eta1p2to1p6 && mereseta_eta1p6to2 &&
0372 mereseta_eta2to2p4) {
0373
0374 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0375
0376
0377 TH1F *resEta1 = mereseta_eta0to0p7->getTH1F();
0378 TH1F *resEta2 = mereseta_eta0p7to1->getTH1F();
0379 TH1F *resEta3 = mereseta_eta1to1p2->getTH1F();
0380 TH1F *resEta4 = mereseta_eta1p2to1p6->getTH1F();
0381 TH1F *resEta5 = mereseta_eta1p6to2->getTH1F();
0382 TH1F *resEta6 = mereseta_eta2to2p4->getTH1F();
0383
0384
0385 MonitorElement *me_res_eta = ibooker.book1D("EtaResolution", "#eta resolution vs |#eta|", eta_binnum, eta_bins);
0386 TH1F *resEta = me_res_eta->getTH1F();
0387 resEta->GetXaxis()->SetTitle("tracking particle |#eta|");
0388 resEta->GetYaxis()->SetTitle("#sigma(#Delta#eta)");
0389 resEta->SetMinimum(0.0);
0390 resEta->SetStats(false);
0391
0392 std::vector<TH1F *> vResEta = {resEta1, resEta2, resEta3, resEta4, resEta5, resEta6};
0393 for (int i = 0; i < 6; i++) {
0394 resEta->SetBinContent(i + 1, vResEta[i]->GetStdDev());
0395 resEta->SetBinError(i + 1, vResEta[i]->GetStdDevError());
0396 }
0397 }
0398 else {
0399 edm::LogWarning("DataNotFound") << "Monitor elements for eta resolution cannot be found!\n";
0400 }
0401
0402 if (meresphi_eta0to0p7 && meresphi_eta0p7to1 && meresphi_eta1to1p2 && meresphi_eta1p2to1p6 && meresphi_eta1p6to2 &&
0403 meresphi_eta2to2p4) {
0404
0405 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0406
0407
0408 TH1F *resPhi1 = meresphi_eta0to0p7->getTH1F();
0409 TH1F *resPhi2 = meresphi_eta0p7to1->getTH1F();
0410 TH1F *resPhi3 = meresphi_eta1to1p2->getTH1F();
0411 TH1F *resPhi4 = meresphi_eta1p2to1p6->getTH1F();
0412 TH1F *resPhi5 = meresphi_eta1p6to2->getTH1F();
0413 TH1F *resPhi6 = meresphi_eta2to2p4->getTH1F();
0414
0415
0416 MonitorElement *me_res_phi = ibooker.book1D("PhiResolution", "#phi resolution vs |#eta|", eta_binnum, eta_bins);
0417 TH1F *resPhi = me_res_phi->getTH1F();
0418 resPhi->GetXaxis()->SetTitle("tracking particle |#eta|");
0419 resPhi->GetYaxis()->SetTitle("#sigma(#Delta#phi)");
0420 resPhi->SetMinimum(0.0);
0421 resPhi->SetStats(false);
0422
0423 std::vector<TH1F *> vResPhi = {resPhi1, resPhi2, resPhi3, resPhi4, resPhi5, resPhi6};
0424 for (int i = 0; i < 6; i++) {
0425 resPhi->SetBinContent(i + 1, vResPhi[i]->GetStdDev());
0426 resPhi->SetBinError(i + 1, vResPhi[i]->GetStdDevError());
0427 }
0428 }
0429 else {
0430 edm::LogWarning("DataNotFound") << "Monitor elements for phi resolution cannot be found!\n";
0431 }
0432
0433 if (meresVtxZ_eta0to0p7 && meresVtxZ_eta0p7to1 && meresVtxZ_eta1to1p2 && meresVtxZ_eta1p2to1p6 &&
0434 meresVtxZ_eta1p6to2 && meresVtxZ_eta2to2p4) {
0435
0436 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0437
0438
0439 TH1F *resVtxZ_1 = meresVtxZ_eta0to0p7->getTH1F();
0440 TH1F *resVtxZ_2 = meresVtxZ_eta0p7to1->getTH1F();
0441 TH1F *resVtxZ_3 = meresVtxZ_eta1to1p2->getTH1F();
0442 TH1F *resVtxZ_4 = meresVtxZ_eta1p2to1p6->getTH1F();
0443 TH1F *resVtxZ_5 = meresVtxZ_eta1p6to2->getTH1F();
0444 TH1F *resVtxZ_6 = meresVtxZ_eta2to2p4->getTH1F();
0445
0446
0447 MonitorElement *me_res_VtxZ = ibooker.book1D("VtxZResolution", "VtxZ resolution vs |#eta|", eta_binnum, eta_bins);
0448 TH1F *resVtxZ = me_res_VtxZ->getTH1F();
0449 resVtxZ->GetXaxis()->SetTitle("tracking particle |#eta|");
0450 resVtxZ->GetYaxis()->SetTitle("#sigma(#DeltaVtxZ) [cm]");
0451 resVtxZ->SetMinimum(0.0);
0452 resVtxZ->SetStats(false);
0453
0454 std::vector<TH1F *> vResVtxZ = {resVtxZ_1, resVtxZ_2, resVtxZ_3, resVtxZ_4, resVtxZ_5, resVtxZ_6};
0455 for (int i = 0; i < 6; i++) {
0456 resVtxZ->SetBinContent(i + 1, vResVtxZ[i]->GetStdDev());
0457 resVtxZ->SetBinError(i + 1, vResVtxZ[i]->GetStdDevError());
0458 }
0459 }
0460 else {
0461 edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ resolution cannot be found!\n";
0462 }
0463
0464 if (meresd0_eta0to0p7 && meresd0_eta0p7to1 && meresd0_eta1to1p2 && meresd0_eta1p2to1p6 && meresd0_eta1p6to2 &&
0465 meresd0_eta2to2p4) {
0466
0467 dbe->setCurrentFolder(topFolderName_ + "/FinalResolution");
0468
0469
0470 TH1F *resd0_1 = meresd0_eta0to0p7->getTH1F();
0471 TH1F *resd0_2 = meresd0_eta0p7to1->getTH1F();
0472 TH1F *resd0_3 = meresd0_eta1to1p2->getTH1F();
0473 TH1F *resd0_4 = meresd0_eta1p2to1p6->getTH1F();
0474 TH1F *resd0_5 = meresd0_eta1p6to2->getTH1F();
0475 TH1F *resd0_6 = meresd0_eta2to2p4->getTH1F();
0476
0477
0478 MonitorElement *me_res_d0 = ibooker.book1D("d0Resolution", "d_{0} resolution vs |#eta|", eta_binnum, eta_bins);
0479 TH1F *resd0 = me_res_d0->getTH1F();
0480 resd0->GetXaxis()->SetTitle("tracking particle |#eta|");
0481 resd0->GetYaxis()->SetTitle("#sigma(#Deltad_{0}) [cm]");
0482 resd0->SetMinimum(0.0);
0483 resd0->SetStats(false);
0484
0485 std::vector<TH1F *> vResD0 = {resd0_1, resd0_2, resd0_3, resd0_4, resd0_5, resd0_6};
0486 for (int i = 0; i < 6; i++) {
0487 resd0->SetBinContent(i + 1, vResD0[i]->GetStdDev());
0488 resd0->SetBinError(i + 1, vResD0[i]->GetStdDevError());
0489 }
0490 }
0491 else {
0492 edm::LogWarning("DataNotFound") << "Monitor elements for d0 resolution cannot be found!\n";
0493 }
0494
0495 }
0496 else {
0497 edm::LogWarning("DataNotFound") << "Cannot find valid DQM back end \n";
0498 }
0499 }
0500
0501 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0502 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0503 void Phase2OTHarvestTrackingParticles::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0504 edm::ParameterSetDescription desc;
0505 desc.add<std::string>("TopFolderName", "TrackerPhase2OTL1TrackV");
0506 descriptions.add("Phase2OTHarvestTrackingParticles", desc);
0507 }
0508 DEFINE_FWK_MODULE(Phase2OTHarvestTrackingParticles);