Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-20 03:45:59

0001 // Package:    Validation/SiTrackerPhase2V
0002 // Class:      Phase2OTHarvestReconstruction
0003 
0004 /**
0005  * This class is part of the Phase 2 Tracker validation framework and performs
0006  * the harvesting step for tracking particle validation. It processes histograms
0007  * created during the earlier validation steps to calculate efficiencies and
0008  * resolutions for stub reconstruction and tracking performance.
0009  *
0010  * Usage:
0011  * To generate histograms from this code, run the test configuration files provided
0012  * in the DQM/SiTrackerPhase2/test directory. The generated histograms can then be
0013  * analyzed or visualized.
0014  */
0015 
0016 // Updated by: Brandi Skipworth, 2025
0017 
0018 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0019 #include "DQMServices/Core/interface/DQMStore.h"
0020 #include "FWCore/Framework/interface/Event.h"
0021 #include "FWCore/Framework/interface/Frameworkfwd.h"
0022 #include "FWCore/Framework/interface/MakerMacros.h"
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0025 #include "FWCore/ServiceRegistry/interface/Service.h"
0026 
0027 class Phase2OTHarvestReconstruction : public DQMEDHarvester {
0028 public:
0029   explicit Phase2OTHarvestReconstruction(const edm::ParameterSet &);
0030   ~Phase2OTHarvestReconstruction() override;
0031   void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override;
0032   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0033 
0034 private:
0035   // ----------member data ---------------------------
0036   std::string topFolderName_;
0037 };
0038 
0039 Phase2OTHarvestReconstruction::Phase2OTHarvestReconstruction(const edm::ParameterSet &iConfig)
0040     : topFolderName_(iConfig.getParameter<std::string>("TopFolderName")) {}
0041 
0042 Phase2OTHarvestReconstruction::~Phase2OTHarvestReconstruction() {}
0043 
0044 // ------------ method called once each job just after ending the event loop
0045 // ------------
0046 void Phase2OTHarvestReconstruction::dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0047   using namespace edm;
0048 
0049   float eta_bins[] = {0.0, 0.7, 1.0, 1.2, 1.6, 2.0, 2.4};
0050   int eta_binnum = 6;
0051 
0052   // Find all monitor elements for histograms
0053   MonitorElement *meN_clus_barrel = igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_if_stub_barrel");
0054   MonitorElement *meD_clus_barrel = igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_barrel");
0055   MonitorElement *meN_clus_zoom_barrel =
0056       igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_if_stub_zoom_barrel");
0057   MonitorElement *meD_clus_zoom_barrel =
0058       igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_zoom_barrel");
0059   MonitorElement *meN_clus_endcaps =
0060       igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_if_stub_endcaps");
0061   MonitorElement *meD_clus_endcaps = igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_endcaps");
0062   MonitorElement *meN_clus_zoom_endcaps =
0063       igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_if_stub_zoom_endcaps");
0064   MonitorElement *meD_clus_zoom_endcaps =
0065       igetter.get(topFolderName_ + "/EfficiencyIngredients/gen_clusters_zoom_endcaps");
0066   MonitorElement *meN_eta = igetter.get(topFolderName_ + "/EfficiencyIngredients/match_tp_eta");
0067   MonitorElement *meD_eta = igetter.get(topFolderName_ + "/EfficiencyIngredients/tp_eta");
0068   MonitorElement *meN_pt = igetter.get(topFolderName_ + "/EfficiencyIngredients/match_tp_pt");
0069   MonitorElement *meD_pt = igetter.get(topFolderName_ + "/EfficiencyIngredients/tp_pt");
0070   MonitorElement *meN_pt_zoom = igetter.get(topFolderName_ + "/EfficiencyIngredients/match_tp_pt_zoom");
0071   MonitorElement *meD_pt_zoom = igetter.get(topFolderName_ + "/EfficiencyIngredients/tp_pt_zoom");
0072   MonitorElement *meN_d0 = igetter.get(topFolderName_ + "/EfficiencyIngredients/match_tp_d0");
0073   MonitorElement *meD_d0 = igetter.get(topFolderName_ + "/EfficiencyIngredients/tp_d0");
0074   MonitorElement *meN_VtxR = igetter.get(topFolderName_ + "/EfficiencyIngredients/match_tp_VtxR");
0075   MonitorElement *meD_VtxR = igetter.get(topFolderName_ + "/EfficiencyIngredients/tp_VtxR");
0076   MonitorElement *meN_VtxZ = igetter.get(topFolderName_ + "/EfficiencyIngredients/match_tp_VtxZ");
0077   MonitorElement *meD_VtxZ = igetter.get(topFolderName_ + "/EfficiencyIngredients/tp_VtxZ");
0078 
0079   MonitorElement *merespt_eta0to0p7_pt2to3 =
0080       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta0to0p7_pt2to3");
0081   MonitorElement *merespt_eta0p7to1_pt2to3 =
0082       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta0p7to1_pt2to3");
0083   MonitorElement *merespt_eta1to1p2_pt2to3 =
0084       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1to1p2_pt2to3");
0085   MonitorElement *merespt_eta1p2to1p6_pt2to3 =
0086       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1p2to1p6_pt2to3");
0087   MonitorElement *merespt_eta1p6to2_pt2to3 =
0088       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1p6to2_pt2to3");
0089   MonitorElement *merespt_eta2to2p4_pt2to3 =
0090       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta2to2p4_pt2to3");
0091   MonitorElement *merespt_eta0to0p7_pt3to8 =
0092       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta0to0p7_pt3to8");
0093   MonitorElement *merespt_eta0p7to1_pt3to8 =
0094       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta0p7to1_pt3to8");
0095   MonitorElement *merespt_eta1to1p2_pt3to8 =
0096       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1to1p2_pt3to8");
0097   MonitorElement *merespt_eta1p2to1p6_pt3to8 =
0098       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1p2to1p6_pt3to8");
0099   MonitorElement *merespt_eta1p6to2_pt3to8 =
0100       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1p6to2_pt3to8");
0101   MonitorElement *merespt_eta2to2p4_pt3to8 =
0102       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta2to2p4_pt3to8");
0103   MonitorElement *merespt_eta0to0p7_pt8toInf =
0104       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta0to0p7_pt8toInf");
0105   MonitorElement *merespt_eta0p7to1_pt8toInf =
0106       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta0p7to1_pt8toInf");
0107   MonitorElement *merespt_eta1to1p2_pt8toInf =
0108       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1to1p2_pt8toInf");
0109   MonitorElement *merespt_eta1p2to1p6_pt8toInf =
0110       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1p2to1p6_pt8toInf");
0111   MonitorElement *merespt_eta1p6to2_pt8toInf =
0112       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta1p6to2_pt8toInf");
0113   MonitorElement *merespt_eta2to2p4_pt8toInf =
0114       igetter.get(topFolderName_ + "/ResolutionIngredients/respt_eta2to2p4_pt8toInf");
0115 
0116   MonitorElement *mereseta_eta0to0p7 = igetter.get(topFolderName_ + "/ResolutionIngredients/reseta_eta0to0p7");
0117   MonitorElement *mereseta_eta0p7to1 = igetter.get(topFolderName_ + "/ResolutionIngredients/reseta_eta0p7to1");
0118   MonitorElement *mereseta_eta1to1p2 = igetter.get(topFolderName_ + "/ResolutionIngredients/reseta_eta1to1p2");
0119   MonitorElement *mereseta_eta1p2to1p6 = igetter.get(topFolderName_ + "/ResolutionIngredients/reseta_eta1p2to1p6");
0120   MonitorElement *mereseta_eta1p6to2 = igetter.get(topFolderName_ + "/ResolutionIngredients/reseta_eta1p6to2");
0121   MonitorElement *mereseta_eta2to2p4 = igetter.get(topFolderName_ + "/ResolutionIngredients/reseta_eta2to2p4");
0122 
0123   MonitorElement *meresphi_eta0to0p7 = igetter.get(topFolderName_ + "/ResolutionIngredients/resphi_eta0to0p7");
0124   MonitorElement *meresphi_eta0p7to1 = igetter.get(topFolderName_ + "/ResolutionIngredients/resphi_eta0p7to1");
0125   MonitorElement *meresphi_eta1to1p2 = igetter.get(topFolderName_ + "/ResolutionIngredients/resphi_eta1to1p2");
0126   MonitorElement *meresphi_eta1p2to1p6 = igetter.get(topFolderName_ + "/ResolutionIngredients/resphi_eta1p2to1p6");
0127   MonitorElement *meresphi_eta1p6to2 = igetter.get(topFolderName_ + "/ResolutionIngredients/resphi_eta1p6to2");
0128   MonitorElement *meresphi_eta2to2p4 = igetter.get(topFolderName_ + "/ResolutionIngredients/resphi_eta2to2p4");
0129 
0130   MonitorElement *meresVtxZ_eta0to0p7 = igetter.get(topFolderName_ + "/ResolutionIngredients/resVtxZ_eta0to0p7");
0131   MonitorElement *meresVtxZ_eta0p7to1 = igetter.get(topFolderName_ + "/ResolutionIngredients/resVtxZ_eta0p7to1");
0132   MonitorElement *meresVtxZ_eta1to1p2 = igetter.get(topFolderName_ + "/ResolutionIngredients/resVtxZ_eta1to1p2");
0133   MonitorElement *meresVtxZ_eta1p2to1p6 = igetter.get(topFolderName_ + "/ResolutionIngredients/resVtxZ_eta1p2to1p6");
0134   MonitorElement *meresVtxZ_eta1p6to2 = igetter.get(topFolderName_ + "/ResolutionIngredients/resVtxZ_eta1p6to2");
0135   MonitorElement *meresVtxZ_eta2to2p4 = igetter.get(topFolderName_ + "/ResolutionIngredients/resVtxZ_eta2to2p4");
0136 
0137   MonitorElement *meresd0_eta0to0p7 = igetter.get(topFolderName_ + "/ResolutionIngredients/resd0_eta0to0p7");
0138   MonitorElement *meresd0_eta0p7to1 = igetter.get(topFolderName_ + "/ResolutionIngredients/resd0_eta0p7to1");
0139   MonitorElement *meresd0_eta1to1p2 = igetter.get(topFolderName_ + "/ResolutionIngredients/resd0_eta1to1p2");
0140   MonitorElement *meresd0_eta1p2to1p6 = igetter.get(topFolderName_ + "/ResolutionIngredients/resd0_eta1p2to1p6");
0141   MonitorElement *meresd0_eta1p6to2 = igetter.get(topFolderName_ + "/ResolutionIngredients/resd0_eta1p6to2");
0142   MonitorElement *meresd0_eta2to2p4 = igetter.get(topFolderName_ + "/ResolutionIngredients/resd0_eta2to2p4");
0143 
0144   if (meN_clus_barrel && meD_clus_barrel) {
0145     // Get the numerator and denominator histograms
0146     TH1F *numerator = meN_clus_barrel->getTH1F();
0147     TH1F *denominator = meD_clus_barrel->getTH1F();
0148     numerator->Sumw2();
0149     denominator->Sumw2();
0150 
0151     // Set the current directory
0152     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0153 
0154     // Book the new histogram to contain the results
0155     MonitorElement *me_effic_clus_barrel = ibooker.book1D("StubEfficiencyBarrel",
0156                                                           "Stub Efficiency Barrel",
0157                                                           numerator->GetNbinsX(),
0158                                                           numerator->GetXaxis()->GetXmin(),
0159                                                           numerator->GetXaxis()->GetXmax());
0160 
0161     // Calculate the efficiency
0162     me_effic_clus_barrel->getTH1F()->Divide(numerator, denominator, 1., 1., "B");
0163     me_effic_clus_barrel->setAxisTitle("tracking particle pT [GeV]");
0164     me_effic_clus_barrel->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0165     me_effic_clus_barrel->getTH1F()->SetMaximum(1.1);
0166     me_effic_clus_barrel->getTH1F()->SetMinimum(0.0);
0167     me_effic_clus_barrel->getTH1F()->SetStats(false);
0168   }  // if ME found
0169   else {
0170     edm::LogWarning("DataNotFound") << "Monitor elements for stub efficiency barrel cannot be found!\n";
0171   }
0172 
0173   if (meN_clus_zoom_barrel && meD_clus_zoom_barrel) {
0174     // Get the numerator and denominator histograms
0175     TH1F *numerator_zoom = meN_clus_zoom_barrel->getTH1F();
0176     TH1F *denominator_zoom = meD_clus_zoom_barrel->getTH1F();
0177     numerator_zoom->Sumw2();
0178     denominator_zoom->Sumw2();
0179 
0180     // Set the current directory
0181     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0182 
0183     // Book the new histogram to contain the results
0184     MonitorElement *me_effic_clus_zoom_barrel = ibooker.book1D("StubEfficiencyZoomBarrel",
0185                                                                "Stub Efficiency Zoom Barrel",
0186                                                                numerator_zoom->GetNbinsX(),
0187                                                                numerator_zoom->GetXaxis()->GetXmin(),
0188                                                                numerator_zoom->GetXaxis()->GetXmax());
0189 
0190     // Calculate the efficiency
0191     me_effic_clus_zoom_barrel->getTH1F()->Divide(numerator_zoom, denominator_zoom, 1., 1., "B");
0192     me_effic_clus_zoom_barrel->setAxisTitle("tracking particle pT [GeV]");
0193     me_effic_clus_zoom_barrel->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0194     me_effic_clus_zoom_barrel->getTH1F()->SetMaximum(1.1);
0195     me_effic_clus_zoom_barrel->getTH1F()->SetMinimum(0.0);
0196     me_effic_clus_zoom_barrel->getTH1F()->SetStats(false);
0197   }  // if ME found
0198   else {
0199     edm::LogWarning("DataNotFound") << "Monitor elements for stub zoom barrel efficiency cannot be found!\n";
0200   }
0201 
0202   if (meN_clus_endcaps && meD_clus_endcaps) {
0203     // Get the numerator and denominator histograms
0204     TH1F *numerator = meN_clus_endcaps->getTH1F();
0205     TH1F *denominator = meD_clus_endcaps->getTH1F();
0206     numerator->Sumw2();
0207     denominator->Sumw2();
0208 
0209     // Set the current directory
0210     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0211 
0212     // Book the new histogram to contain the results
0213     MonitorElement *me_effic_clus_endcaps = ibooker.book1D("StubEfficiencyEndcaps",
0214                                                            "Stub Efficiency Endcaps",
0215                                                            numerator->GetNbinsX(),
0216                                                            numerator->GetXaxis()->GetXmin(),
0217                                                            numerator->GetXaxis()->GetXmax());
0218 
0219     // Calculate the efficiency
0220     me_effic_clus_endcaps->getTH1F()->Divide(numerator, denominator, 1., 1., "B");
0221     me_effic_clus_endcaps->setAxisTitle("tracking particle pT [GeV]");
0222     me_effic_clus_endcaps->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0223     me_effic_clus_endcaps->getTH1F()->SetMaximum(1.1);
0224     me_effic_clus_endcaps->getTH1F()->SetMinimum(0.0);
0225     me_effic_clus_endcaps->getTH1F()->SetStats(false);
0226   }  // if ME found
0227   else {
0228     edm::LogWarning("DataNotFound") << "Monitor elements for stub efficiency endcaps cannot be found!\n";
0229   }
0230 
0231   if (meN_clus_zoom_endcaps && meD_clus_zoom_endcaps) {
0232     // Get the numerator and denominator histograms
0233     TH1F *numerator_zoom = meN_clus_zoom_endcaps->getTH1F();
0234     TH1F *denominator_zoom = meD_clus_zoom_endcaps->getTH1F();
0235     numerator_zoom->Sumw2();
0236     denominator_zoom->Sumw2();
0237 
0238     // Set the current directory
0239     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0240 
0241     // Book the new histogram to contain the results
0242     MonitorElement *me_effic_clus_zoom_endcaps = ibooker.book1D("StubEfficiencyZoomEndcaps",
0243                                                                 "Stub Efficiency Zoom Endcaps",
0244                                                                 numerator_zoom->GetNbinsX(),
0245                                                                 numerator_zoom->GetXaxis()->GetXmin(),
0246                                                                 numerator_zoom->GetXaxis()->GetXmax());
0247 
0248     // Calculate the efficiency
0249     me_effic_clus_zoom_endcaps->getTH1F()->Divide(numerator_zoom, denominator_zoom, 1., 1., "B");
0250     me_effic_clus_zoom_endcaps->setAxisTitle("tracking particle pT [GeV]");
0251     me_effic_clus_zoom_endcaps->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0252     me_effic_clus_zoom_endcaps->getTH1F()->SetMaximum(1.1);
0253     me_effic_clus_zoom_endcaps->getTH1F()->SetMinimum(0.0);
0254     me_effic_clus_zoom_endcaps->getTH1F()->SetStats(false);
0255   }  // if ME found
0256   else {
0257     edm::LogWarning("DataNotFound") << "Monitor elements for stub zoom endcaps efficiency cannot be found!\n";
0258   }
0259 
0260   if (meN_eta && meD_eta) {
0261     // Get the numerator and denominator histograms
0262     TH1F *numerator = meN_eta->getTH1F();
0263     TH1F *denominator = meD_eta->getTH1F();
0264     numerator->Sumw2();
0265     denominator->Sumw2();
0266 
0267     // Set the current directory
0268     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0269 
0270     // Book the new histogram to contain the results
0271     MonitorElement *me_effic_eta = ibooker.book1D("EtaEfficiency",
0272                                                   "#eta efficiency",
0273                                                   numerator->GetNbinsX(),
0274                                                   numerator->GetXaxis()->GetXmin(),
0275                                                   numerator->GetXaxis()->GetXmax());
0276 
0277     // Calculate the efficiency
0278     me_effic_eta->getTH1F()->Divide(numerator, denominator, 1., 1., "B");
0279     me_effic_eta->setAxisTitle("tracking particle #eta");
0280     me_effic_eta->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0281     me_effic_eta->getTH1F()->SetMaximum(1.0);
0282     me_effic_eta->getTH1F()->SetMinimum(0.0);
0283     me_effic_eta->getTH1F()->SetStats(false);
0284   }  // if ME found
0285   else {
0286     edm::LogWarning("DataNotFound") << "Monitor elements for eta efficiency cannot be found!\n";
0287   }
0288 
0289   if (meN_pt && meD_pt) {
0290     // Get the numerator and denominator histograms
0291     TH1F *numerator2 = meN_pt->getTH1F();
0292     numerator2->Sumw2();
0293     TH1F *denominator2 = meD_pt->getTH1F();
0294     denominator2->Sumw2();
0295 
0296     // Set the current directory
0297     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0298 
0299     // Book the new histogram to contain the results
0300     MonitorElement *me_effic_pt = ibooker.book1D("PtEfficiency",
0301                                                  "p_{T} efficiency",
0302                                                  numerator2->GetNbinsX(),
0303                                                  numerator2->GetXaxis()->GetXmin(),
0304                                                  numerator2->GetXaxis()->GetXmax());
0305 
0306     // Calculate the efficiency
0307     me_effic_pt->getTH1F()->Divide(numerator2, denominator2, 1., 1., "B");
0308     me_effic_pt->setAxisTitle("Tracking particle p_{T} [GeV]");
0309     me_effic_pt->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0310     me_effic_pt->getTH1F()->SetMaximum(1.0);
0311     me_effic_pt->getTH1F()->SetMinimum(0.0);
0312     me_effic_pt->getTH1F()->SetStats(false);
0313   }  // if ME found
0314   else {
0315     edm::LogWarning("DataNotFound") << "Monitor elements for pT efficiency cannot be found!\n";
0316   }
0317 
0318   if (meN_pt_zoom && meD_pt_zoom) {
0319     // Get the numerator and denominator histograms
0320     TH1F *numerator2_zoom = meN_pt_zoom->getTH1F();
0321     numerator2_zoom->Sumw2();
0322     TH1F *denominator2_zoom = meD_pt_zoom->getTH1F();
0323     denominator2_zoom->Sumw2();
0324 
0325     // Set the current directory
0326     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0327 
0328     // Book the new histogram to contain the results
0329     MonitorElement *me_effic_pt_zoom = ibooker.book1D("PtEfficiency_zoom",
0330                                                       "p_{T} efficiency",
0331                                                       numerator2_zoom->GetNbinsX(),
0332                                                       numerator2_zoom->GetXaxis()->GetXmin(),
0333                                                       numerator2_zoom->GetXaxis()->GetXmax());
0334 
0335     // Calculate the efficiency
0336     me_effic_pt_zoom->getTH1F()->Divide(numerator2_zoom, denominator2_zoom, 1., 1., "B");
0337     me_effic_pt_zoom->setAxisTitle("Tracking particle p_{T} [GeV]");
0338     me_effic_pt_zoom->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0339     me_effic_pt_zoom->getTH1F()->SetMaximum(1.0);
0340     me_effic_pt_zoom->getTH1F()->SetMinimum(0.0);
0341     me_effic_pt_zoom->getTH1F()->SetStats(false);
0342   }  // if ME found
0343   else {
0344     edm::LogWarning("DataNotFound") << "Monitor elements for zoom pT efficiency cannot be found!\n";
0345   }
0346 
0347   if (meN_d0 && meD_d0) {
0348     // Get the numerator and denominator histograms
0349     TH1F *numerator5 = meN_d0->getTH1F();
0350     numerator5->Sumw2();
0351     TH1F *denominator5 = meD_d0->getTH1F();
0352     denominator5->Sumw2();
0353 
0354     // Set the current directory
0355     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0356 
0357     // Book the new histogram to contain the results
0358     MonitorElement *me_effic_d0 = ibooker.book1D("d0Efficiency",
0359                                                  "d_{0} efficiency",
0360                                                  numerator5->GetNbinsX(),
0361                                                  numerator5->GetXaxis()->GetXmin(),
0362                                                  numerator5->GetXaxis()->GetXmax());
0363 
0364     // Calculate the efficiency
0365     me_effic_d0->getTH1F()->Divide(numerator5, denominator5, 1., 1., "B");
0366     me_effic_d0->setAxisTitle("Tracking particle d_{0} [cm]");
0367     me_effic_d0->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0368     me_effic_d0->getTH1F()->SetMaximum(1.0);
0369     me_effic_d0->getTH1F()->SetMinimum(0.0);
0370     me_effic_d0->getTH1F()->SetStats(false);
0371   }  // if ME found
0372   else {
0373     edm::LogWarning("DataNotFound") << "Monitor elements for d0 efficiency cannot be found!\n";
0374   }
0375 
0376   if (meN_VtxR && meD_VtxR) {
0377     // Get the numerator and denominator histograms
0378     TH1F *numerator6 = meN_VtxR->getTH1F();
0379     numerator6->Sumw2();
0380     TH1F *denominator6 = meD_VtxR->getTH1F();
0381     denominator6->Sumw2();
0382 
0383     // Set the current directory
0384     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0385 
0386     // Book the new histogram to contain the results
0387     MonitorElement *me_effic_VtxR = ibooker.book1D("VtxREfficiency",
0388                                                    "Vtx R efficiency",
0389                                                    numerator6->GetNbinsX(),
0390                                                    numerator6->GetXaxis()->GetXmin(),
0391                                                    numerator6->GetXaxis()->GetXmax());
0392 
0393     // Calculate the efficiency
0394     me_effic_VtxR->getTH1F()->Divide(numerator6, denominator6, 1., 1., "B");
0395     me_effic_VtxR->setAxisTitle("Tracking particle VtxR [cm]");
0396     me_effic_VtxR->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0397     me_effic_VtxR->getTH1F()->SetMaximum(1.0);
0398     me_effic_VtxR->getTH1F()->SetMinimum(0.0);
0399     me_effic_VtxR->getTH1F()->SetStats(false);
0400   }  // if ME found
0401   else {
0402     edm::LogWarning("DataNotFound") << "Monitor elements for VtxR efficiency cannot be found!\n";
0403   }
0404 
0405   if (meN_VtxZ && meD_VtxZ) {
0406     // Get the numerator and denominator histograms
0407     TH1F *numerator7 = meN_VtxZ->getTH1F();
0408     numerator7->Sumw2();
0409     TH1F *denominator7 = meD_VtxZ->getTH1F();
0410     denominator7->Sumw2();
0411 
0412     // Set the current directory
0413     igetter.setCurrentFolder(topFolderName_ + "/FinalEfficiency");
0414 
0415     // Book the new histogram to contain the results
0416     MonitorElement *me_effic_VtxZ = ibooker.book1D("VtxZEfficiency",
0417                                                    "Vtx Z efficiency",
0418                                                    numerator7->GetNbinsX(),
0419                                                    numerator7->GetXaxis()->GetXmin(),
0420                                                    numerator7->GetXaxis()->GetXmax());
0421 
0422     // Calculate the efficiency
0423     me_effic_VtxZ->getTH1F()->Divide(numerator7, denominator7, 1., 1., "B");
0424     me_effic_VtxZ->setAxisTitle("Tracking particle VtxZ [cm]");
0425     me_effic_VtxZ->getTH1F()->GetYaxis()->SetTitle("Efficiency");
0426     me_effic_VtxZ->getTH1F()->SetMaximum(1.0);
0427     me_effic_VtxZ->getTH1F()->SetMinimum(0.0);
0428     me_effic_VtxZ->getTH1F()->SetStats(false);
0429   }  // if ME found
0430   else {
0431     edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ efficiency cannot be found!\n";
0432   }
0433 
0434   if (merespt_eta0to0p7_pt2to3 && merespt_eta0p7to1_pt2to3 && merespt_eta1to1p2_pt2to3 && merespt_eta1p2to1p6_pt2to3 &&
0435       merespt_eta1p6to2_pt2to3 && merespt_eta2to2p4_pt2to3) {
0436     // Set the current directoy
0437     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0438 
0439     // Grab the histograms
0440     TH1F *resPt1a = merespt_eta0to0p7_pt2to3->getTH1F();
0441     TH1F *resPt2a = merespt_eta0p7to1_pt2to3->getTH1F();
0442     TH1F *resPt3a = merespt_eta1to1p2_pt2to3->getTH1F();
0443     TH1F *resPt4a = merespt_eta1p2to1p6_pt2to3->getTH1F();
0444     TH1F *resPt5a = merespt_eta1p6to2_pt2to3->getTH1F();
0445     TH1F *resPt6a = merespt_eta2to2p4_pt2to3->getTH1F();
0446 
0447     // Book the new histogram to contain the results
0448     MonitorElement *me_res_pt1 =
0449         ibooker.book1D("pTResVsEta_2-3", "p_{T} resolution vs |#eta|, for p_{T}: 2-3 GeV", eta_binnum, eta_bins);
0450     TH1F *resPt1 = me_res_pt1->getTH1F();
0451     resPt1->GetXaxis()->SetTitle("tracking particle |#eta|");
0452     resPt1->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
0453     resPt1->SetMinimum(0.0);
0454     resPt1->SetStats(false);
0455 
0456     std::vector<TH1F *> vResPt1 = {resPt1a, resPt2a, resPt3a, resPt4a, resPt5a, resPt6a};
0457     for (int i = 0; i < 6; i++) {
0458       resPt1->SetBinContent(i + 1, vResPt1[i]->GetStdDev());
0459       resPt1->SetBinError(i + 1, vResPt1[i]->GetStdDevError());
0460     }
0461   }  // if ME found
0462   else {
0463     edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (2-3) cannot be found!\n";
0464   }
0465 
0466   if (merespt_eta0to0p7_pt3to8 && merespt_eta0p7to1_pt3to8 && merespt_eta1to1p2_pt3to8 && merespt_eta1p2to1p6_pt3to8 &&
0467       merespt_eta1p6to2_pt3to8 && merespt_eta2to2p4_pt3to8) {
0468     // Set the current directoy
0469     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0470 
0471     // Grab the histograms
0472     TH1F *resPt1b = merespt_eta0to0p7_pt3to8->getTH1F();
0473     TH1F *resPt2b = merespt_eta0p7to1_pt3to8->getTH1F();
0474     TH1F *resPt3b = merespt_eta1to1p2_pt3to8->getTH1F();
0475     TH1F *resPt4b = merespt_eta1p2to1p6_pt3to8->getTH1F();
0476     TH1F *resPt5b = merespt_eta1p6to2_pt3to8->getTH1F();
0477     TH1F *resPt6b = merespt_eta2to2p4_pt3to8->getTH1F();
0478 
0479     // Book the new histogram to contain the results
0480     MonitorElement *me_res_pt2 =
0481         ibooker.book1D("pTResVsEta_3-8", "p_{T} resolution vs |#eta|, for p_{T}: 3-8 GeV", eta_binnum, eta_bins);
0482     TH1F *resPt2 = me_res_pt2->getTH1F();
0483     resPt2->GetXaxis()->SetTitle("tracking particle |#eta|");
0484     resPt2->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
0485     resPt2->SetMinimum(0.0);
0486     resPt2->SetStats(false);
0487 
0488     std::vector<TH1F *> vResPt2 = {resPt1b, resPt2b, resPt3b, resPt4b, resPt5b, resPt6b};
0489     for (int i = 0; i < 6; i++) {
0490       resPt2->SetBinContent(i + 1, vResPt2[i]->GetStdDev());
0491       resPt2->SetBinError(i + 1, vResPt2[i]->GetStdDevError());
0492     }
0493   }  // if ME found
0494   else {
0495     edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (3-8) cannot be found!\n";
0496   }
0497 
0498   if (merespt_eta0to0p7_pt8toInf && merespt_eta0p7to1_pt8toInf && merespt_eta1to1p2_pt8toInf &&
0499       merespt_eta1p2to1p6_pt8toInf && merespt_eta1p6to2_pt8toInf && merespt_eta2to2p4_pt8toInf) {
0500     // Set the current directoy
0501     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0502 
0503     // Grab the histograms
0504     TH1F *resPt1c = merespt_eta0to0p7_pt8toInf->getTH1F();
0505     TH1F *resPt2c = merespt_eta0p7to1_pt8toInf->getTH1F();
0506     TH1F *resPt3c = merespt_eta1to1p2_pt8toInf->getTH1F();
0507     TH1F *resPt4c = merespt_eta1p2to1p6_pt8toInf->getTH1F();
0508     TH1F *resPt5c = merespt_eta1p6to2_pt8toInf->getTH1F();
0509     TH1F *resPt6c = merespt_eta2to2p4_pt8toInf->getTH1F();
0510 
0511     // Book the new histogram to contain the results
0512     MonitorElement *me_res_pt3 =
0513         ibooker.book1D("pTResVsEta_8-inf", "p_{T} resolution vs |#eta|, for p_{T}: >8 GeV", eta_binnum, eta_bins);
0514     TH1F *resPt3 = me_res_pt3->getTH1F();
0515     resPt3->GetXaxis()->SetTitle("tracking particle |#eta|");
0516     resPt3->GetYaxis()->SetTitle("#sigma(#Deltap_{T}/p_{T})");
0517     resPt3->SetMinimum(0.0);
0518     resPt3->SetStats(false);
0519 
0520     std::vector<TH1F *> vResPt3 = {resPt1c, resPt2c, resPt3c, resPt4c, resPt5c, resPt6c};
0521     for (int i = 0; i < 6; i++) {
0522       resPt3->SetBinContent(i + 1, vResPt3[i]->GetStdDev());
0523       resPt3->SetBinError(i + 1, vResPt3[i]->GetStdDevError());
0524     }
0525   }  // if ME found
0526   else {
0527     edm::LogWarning("DataNotFound") << "Monitor elements for pT resolution (8-inf) cannot be found!\n";
0528   }
0529 
0530   if (mereseta_eta0to0p7 && mereseta_eta0p7to1 && mereseta_eta1to1p2 && mereseta_eta1p2to1p6 && mereseta_eta1p6to2 &&
0531       mereseta_eta2to2p4) {
0532     // Set the current directoy
0533     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0534 
0535     // Grab the histograms
0536     TH1F *resEta1 = mereseta_eta0to0p7->getTH1F();
0537     TH1F *resEta2 = mereseta_eta0p7to1->getTH1F();
0538     TH1F *resEta3 = mereseta_eta1to1p2->getTH1F();
0539     TH1F *resEta4 = mereseta_eta1p2to1p6->getTH1F();
0540     TH1F *resEta5 = mereseta_eta1p6to2->getTH1F();
0541     TH1F *resEta6 = mereseta_eta2to2p4->getTH1F();
0542 
0543     // Book the new histogram to contain the results
0544     MonitorElement *me_res_eta = ibooker.book1D("EtaResolution", "#eta resolution vs |#eta|", eta_binnum, eta_bins);
0545     TH1F *resEta = me_res_eta->getTH1F();
0546     resEta->GetXaxis()->SetTitle("tracking particle |#eta|");
0547     resEta->GetYaxis()->SetTitle("#sigma(#Delta#eta)");
0548     resEta->SetMinimum(0.0);
0549     resEta->SetStats(false);
0550 
0551     std::vector<TH1F *> vResEta = {resEta1, resEta2, resEta3, resEta4, resEta5, resEta6};
0552     for (int i = 0; i < 6; i++) {
0553       resEta->SetBinContent(i + 1, vResEta[i]->GetStdDev());
0554       resEta->SetBinError(i + 1, vResEta[i]->GetStdDevError());
0555     }
0556   }  // if ME found
0557   else {
0558     edm::LogWarning("DataNotFound") << "Monitor elements for eta resolution cannot be found!\n";
0559   }
0560 
0561   if (meresphi_eta0to0p7 && meresphi_eta0p7to1 && meresphi_eta1to1p2 && meresphi_eta1p2to1p6 && meresphi_eta1p6to2 &&
0562       meresphi_eta2to2p4) {
0563     // Set the current directoy
0564     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0565 
0566     // Grab the histograms
0567     TH1F *resPhi1 = meresphi_eta0to0p7->getTH1F();
0568     TH1F *resPhi2 = meresphi_eta0p7to1->getTH1F();
0569     TH1F *resPhi3 = meresphi_eta1to1p2->getTH1F();
0570     TH1F *resPhi4 = meresphi_eta1p2to1p6->getTH1F();
0571     TH1F *resPhi5 = meresphi_eta1p6to2->getTH1F();
0572     TH1F *resPhi6 = meresphi_eta2to2p4->getTH1F();
0573 
0574     // Book the new histogram to contain the results
0575     MonitorElement *me_res_phi = ibooker.book1D("PhiResolution", "#phi resolution vs |#eta|", eta_binnum, eta_bins);
0576     TH1F *resPhi = me_res_phi->getTH1F();
0577     resPhi->GetXaxis()->SetTitle("tracking particle |#eta|");
0578     resPhi->GetYaxis()->SetTitle("#sigma(#Delta#phi)");
0579     resPhi->SetMinimum(0.0);
0580     resPhi->SetStats(false);
0581 
0582     std::vector<TH1F *> vResPhi = {resPhi1, resPhi2, resPhi3, resPhi4, resPhi5, resPhi6};
0583     for (int i = 0; i < 6; i++) {
0584       resPhi->SetBinContent(i + 1, vResPhi[i]->GetStdDev());
0585       resPhi->SetBinError(i + 1, vResPhi[i]->GetStdDevError());
0586     }
0587   }  // if ME found
0588   else {
0589     edm::LogWarning("DataNotFound") << "Monitor elements for phi resolution cannot be found!\n";
0590   }
0591 
0592   if (meresVtxZ_eta0to0p7 && meresVtxZ_eta0p7to1 && meresVtxZ_eta1to1p2 && meresVtxZ_eta1p2to1p6 &&
0593       meresVtxZ_eta1p6to2 && meresVtxZ_eta2to2p4) {
0594     // Set the current directoy
0595     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0596 
0597     // Grab the histograms
0598     TH1F *resVtxZ_1 = meresVtxZ_eta0to0p7->getTH1F();
0599     TH1F *resVtxZ_2 = meresVtxZ_eta0p7to1->getTH1F();
0600     TH1F *resVtxZ_3 = meresVtxZ_eta1to1p2->getTH1F();
0601     TH1F *resVtxZ_4 = meresVtxZ_eta1p2to1p6->getTH1F();
0602     TH1F *resVtxZ_5 = meresVtxZ_eta1p6to2->getTH1F();
0603     TH1F *resVtxZ_6 = meresVtxZ_eta2to2p4->getTH1F();
0604 
0605     // Book the new histogram to contain the results
0606     MonitorElement *me_res_VtxZ = ibooker.book1D("VtxZResolution", "VtxZ resolution vs |#eta|", eta_binnum, eta_bins);
0607     TH1F *resVtxZ = me_res_VtxZ->getTH1F();
0608     resVtxZ->GetXaxis()->SetTitle("tracking particle |#eta|");
0609     resVtxZ->GetYaxis()->SetTitle("#sigma(#DeltaVtxZ) [cm]");
0610     resVtxZ->SetMinimum(0.0);
0611     resVtxZ->SetStats(false);
0612 
0613     std::vector<TH1F *> vResVtxZ = {resVtxZ_1, resVtxZ_2, resVtxZ_3, resVtxZ_4, resVtxZ_5, resVtxZ_6};
0614     for (int i = 0; i < 6; i++) {
0615       resVtxZ->SetBinContent(i + 1, vResVtxZ[i]->GetStdDev());
0616       resVtxZ->SetBinError(i + 1, vResVtxZ[i]->GetStdDevError());
0617     }
0618   }  // if ME found
0619   else {
0620     edm::LogWarning("DataNotFound") << "Monitor elements for VtxZ resolution cannot be found!\n";
0621   }
0622 
0623   if (meresd0_eta0to0p7 && meresd0_eta0p7to1 && meresd0_eta1to1p2 && meresd0_eta1p2to1p6 && meresd0_eta1p6to2 &&
0624       meresd0_eta2to2p4) {
0625     // Set the current directoy
0626     igetter.setCurrentFolder(topFolderName_ + "/FinalResolution");
0627 
0628     // Grab the histograms
0629     TH1F *resd0_1 = meresd0_eta0to0p7->getTH1F();
0630     TH1F *resd0_2 = meresd0_eta0p7to1->getTH1F();
0631     TH1F *resd0_3 = meresd0_eta1to1p2->getTH1F();
0632     TH1F *resd0_4 = meresd0_eta1p2to1p6->getTH1F();
0633     TH1F *resd0_5 = meresd0_eta1p6to2->getTH1F();
0634     TH1F *resd0_6 = meresd0_eta2to2p4->getTH1F();
0635 
0636     // Book the new histogram to contain the results
0637     MonitorElement *me_res_d0 = ibooker.book1D("d0Resolution", "d_{0} resolution vs |#eta|", eta_binnum, eta_bins);
0638     TH1F *resd0 = me_res_d0->getTH1F();
0639     resd0->GetXaxis()->SetTitle("tracking particle |#eta|");
0640     resd0->GetYaxis()->SetTitle("#sigma(#Deltad_{0}) [cm]");
0641     resd0->SetMinimum(0.0);
0642     resd0->SetStats(false);
0643 
0644     std::vector<TH1F *> vResD0 = {resd0_1, resd0_2, resd0_3, resd0_4, resd0_5, resd0_6};
0645     for (int i = 0; i < 6; i++) {
0646       resd0->SetBinContent(i + 1, vResD0[i]->GetStdDev());
0647       resd0->SetBinError(i + 1, vResD0[i]->GetStdDevError());
0648     }
0649   }  // if ME found
0650   else {
0651     edm::LogWarning("DataNotFound") << "Monitor elements for d0 resolution cannot be found!\n";
0652   }
0653 }  // end dqmEndJob
0654 
0655 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0656 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0657 void Phase2OTHarvestReconstruction::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0658   edm::ParameterSetDescription desc;
0659   desc.add<std::string>("TopFolderName", "TrackerPhase2OTL1TrackV");
0660   descriptions.add("Phase2OTHarvestReconstruction", desc);
0661 }
0662 DEFINE_FWK_MODULE(Phase2OTHarvestReconstruction);