Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:48

0001 #include "HLTriggerOffline/Btag/interface/HLTVertexPerformanceAnalyzer.h"
0002 
0003 using namespace edm;
0004 using namespace reco;
0005 
0006 HLTVertexPerformanceAnalyzer::HLTVertexPerformanceAnalyzer(const edm::ParameterSet &iConfig) {
0007   mainFolder_ = iConfig.getParameter<std::string>("mainFolder");
0008   hlTriggerResults_ = consumes<TriggerResults>(iConfig.getParameter<InputTag>("TriggerResults"));
0009   VertexCollection_ =
0010       edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag>>("Vertex"),
0011                             [this](edm::InputTag const &tag) { return mayConsume<reco::VertexCollection>(tag); });
0012   hltPathNames_ = iConfig.getParameter<std::vector<std::string>>("HLTPathNames");
0013   simVertexCollection_ = consumes<std::vector<SimVertex>>(iConfig.getParameter<edm::InputTag>("SimVertexCollection"));
0014 
0015   EDConsumerBase::labelsForToken(hlTriggerResults_, label);
0016   hlTriggerResults_Label = label.module;
0017 
0018   for (unsigned int i = 0; i < VertexCollection_.size(); i++) {
0019     EDConsumerBase::labelsForToken(VertexCollection_[i], label);
0020     VertexCollection_Label.push_back(label.module);
0021   }
0022 }
0023 
0024 HLTVertexPerformanceAnalyzer::~HLTVertexPerformanceAnalyzer() {}
0025 
0026 void HLTVertexPerformanceAnalyzer::dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) {
0027   triggerConfChanged_ = true;
0028   EDConsumerBase::labelsForToken(hlTriggerResults_, label);
0029 
0030   hltConfigProvider_.init(iRun, iSetup, label.process, triggerConfChanged_);
0031   const std::vector<std::string> &allHltPathNames = hltConfigProvider_.triggerNames();
0032 
0033   // fill hltPathIndexs_ with the trigger number of each hltPathNames_
0034   for (size_t trgs = 0; trgs < hltPathNames_.size(); trgs++) {
0035     unsigned int found = 1;
0036     int it_mem = -1;
0037     for (size_t it = 0; it < allHltPathNames.size(); ++it) {
0038       found = allHltPathNames.at(it).find(hltPathNames_[trgs]);
0039       if (found == 0) {
0040         it_mem = (int)it;
0041       }
0042     }  // for allallHltPathNames
0043     hltPathIndexs_.push_back(it_mem);
0044   }  // for hltPathNames_
0045 
0046   // fill _isfoundHLTs for each hltPathNames_
0047   for (size_t trgs = 0; trgs < hltPathNames_.size(); trgs++) {
0048     if (hltPathIndexs_[trgs] < 0) {
0049       _isfoundHLTs.push_back(false);
0050     } else {
0051       _isfoundHLTs.push_back(true);
0052     }
0053   }
0054 }
0055 
0056 void HLTVertexPerformanceAnalyzer::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0057   bool trigRes = false;
0058   using namespace edm;
0059 
0060   // get triggerResults
0061   Handle<TriggerResults> TriggerResulsHandler;
0062   if (hlTriggerResults_Label.empty() || hlTriggerResults_Label == "NULL") {
0063     edm::LogInfo("NoTriggerResults") << "TriggerResults ==> Empty";
0064     return;
0065   }
0066   iEvent.getByToken(hlTriggerResults_, TriggerResulsHandler);
0067   if (TriggerResulsHandler.isValid())
0068     trigRes = true;
0069   if (!trigRes) {
0070     edm::LogInfo("NoTriggerResults") << "TriggerResults ==> not readable";
0071     return;
0072   }
0073   const TriggerResults &triggerResults = *(TriggerResulsHandler.product());
0074 
0075   // get simVertex
0076   float simPV = 0;
0077 
0078   Handle<std::vector<SimVertex>> simVertexCollection;
0079   iEvent.getByToken(simVertexCollection_, simVertexCollection);
0080   if (!simVertexCollection.isValid()) {
0081     edm::LogInfo("NoSimVertex") << "SimVertex collection is invalid";
0082     return;
0083   }
0084   if (simVertexCollection->empty()) {
0085     edm::LogInfo("NoSimVertex") << "SimVertex collection is empty";
0086     return;
0087   }
0088   const SimVertex simPVh = *(simVertexCollection->begin());
0089   simPV = simPVh.position().z();
0090 
0091   // fill the DQM plot
0092   Handle<VertexCollection> VertexHandler;
0093   for (unsigned int ind = 0; ind < hltPathNames_.size(); ind++) {
0094     for (unsigned int coll = 0; coll < VertexCollection_.size(); coll++) {
0095       bool VertexOK = false;
0096       if (!_isfoundHLTs[ind])
0097         continue;  // if the hltPath is not in the event, skip the event
0098       if (!triggerResults.accept(hltPathIndexs_[ind]))
0099         continue;  // if the hltPath was not accepted skip the event
0100 
0101       // get the recoVertex
0102       if (!VertexCollection_Label.at(coll).empty() && VertexCollection_Label.at(coll) != "NULL") {
0103         iEvent.getByToken(VertexCollection_.at(coll), VertexHandler);
0104         if (VertexHandler.isValid() > 0)
0105           VertexOK = true;
0106       }
0107 
0108       if (VertexOK) {
0109         // calculate the variable (RecoVertex - SimVertex)
0110         float value = VertexHandler->begin()->z() - simPV;
0111 
0112         // if value is over/under flow, assign the extreme value
0113         float maxValue = H1_.at(ind)["Vertex_" + VertexCollection_Label.at(coll)]->getTH1F()->GetXaxis()->GetXmax();
0114         if (value > maxValue)
0115           value = maxValue - 0.0001;
0116         if (value < -maxValue)
0117           value = -maxValue + 0.0001;
0118         // fill the histo
0119         H1_.at(ind)["Vertex_" + VertexCollection_Label.at(coll)]->Fill(value);
0120       }
0121     }  // for on VertexCollection_
0122   }    // for on hltPathNames_
0123 }
0124 
0125 void HLTVertexPerformanceAnalyzer::bookHistograms(DQMStore::IBooker &ibooker,
0126                                                   edm::Run const &iRun,
0127                                                   edm::EventSetup const &iSetup) {
0128   // book the DQM plots
0129   using namespace std;
0130   std::string dqmFolder;
0131   for (unsigned int ind = 0; ind < hltPathNames_.size(); ind++) {
0132     dqmFolder = Form("%s/Vertex/%s", mainFolder_.c_str(), hltPathNames_[ind].c_str());
0133     H1_.push_back(std::map<std::string, MonitorElement *>());
0134     ibooker.setCurrentFolder(dqmFolder);
0135     for (unsigned int coll = 0; coll < VertexCollection_.size(); coll++) {
0136       float maxValue = 0.02;
0137       if (VertexCollection_Label.at(coll) == ("hltFastPrimaryVertex"))
0138         maxValue = 2.;  // for the hltFastPrimaryVertex use a larger scale (res ~ 1 cm)
0139       float vertexU = maxValue;
0140       float vertexL = -maxValue;
0141       int vertexBins = 100;
0142       if (!VertexCollection_Label.at(coll).empty() && VertexCollection_Label.at(coll) != "NULL") {
0143         H1_.back()["Vertex_" + VertexCollection_Label.at(coll)] =
0144             ibooker.book1D("Vertex_" + VertexCollection_Label.at(coll),
0145                            VertexCollection_Label.at(coll).c_str(),
0146                            vertexBins,
0147                            vertexL,
0148                            vertexU);
0149         H1_.back()["Vertex_" + VertexCollection_Label.at(coll)]->setAxisTitle("vertex error (cm)", 1);
0150       }
0151     }
0152   }
0153 }
0154 
0155 DEFINE_FWK_MODULE(HLTVertexPerformanceAnalyzer);