Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-13 02:58:34

0001 #include "Validation/HLTrigger/interface/HLTGenValHistCollFilter.h"
0002 
0003 // constructor
0004 HLTGenValHistCollFilter::HLTGenValHistCollFilter(edm::ParameterSet filterCollConfig) {
0005   objType_ = filterCollConfig.getParameter<std::string>("objType");
0006   filter_ = filterCollConfig.getParameter<std::string>("filterName");
0007   path_ = filterCollConfig.getParameter<std::string>("pathName");
0008   hltProcessName_ = filterCollConfig.getParameter<std::string>("hltProcessName");
0009   dR2limit_ = filterCollConfig.getParameter<double>("dR2limit");
0010   separator_ = "__";
0011 }
0012 
0013 edm::ParameterSetDescription HLTGenValHistCollFilter::makePSetDescription() {
0014   edm::ParameterSetDescription desc;
0015   desc.add<std::string>("objType", "");
0016   desc.add<std::string>("hltProcessName", "HLT");
0017   desc.add<double>("dR2limit", 0.1);
0018   desc.add<std::string>("filterName", "");
0019   desc.add<std::string>("pathName", "");
0020   return desc;
0021 }
0022 
0023 // general hist booking function, receiving configurations for 1D and 2D hists and calling the respective functions
0024 void HLTGenValHistCollFilter::bookHists(DQMStore::IBooker& iBooker,
0025                                         const std::vector<edm::ParameterSet>& histConfigs,
0026                                         const std::vector<edm::ParameterSet>& histConfigs2D) {
0027   for (const auto& histConfig : histConfigs)
0028     book1D(iBooker, histConfig);
0029   for (const auto& histConfig : histConfigs2D)
0030     book2D(iBooker, histConfig);
0031 }
0032 
0033 // histogram filling routine
0034 void HLTGenValHistCollFilter::fillHists(const HLTGenValObject& obj, edm::Handle<trigger::TriggerEvent>& triggerEvent) {
0035   // this handles the "before" step, denoted by a "dummy" filter called "beforeAnyFilter"
0036   // the histogram is filled without any additional requirements for all objects
0037   if (filter_ == "beforeAnyFilter") {
0038     for (auto& hist : hists_)
0039       hist->fill(obj);
0040   } else {
0041     // main filling code
0042 
0043     // get filter object from filter name
0044     // filters that are ignored start with "-" in the path but are not that way in trigger results
0045     const std::string& filterName = !filter_.empty() && filter_[0] == '-' ? filter_.substr(1) : filter_;
0046     edm::InputTag filterTag(filterName, "", hltProcessName_);
0047     size_t filterIndex = triggerEvent->filterIndex(filterTag);
0048     // get trigger objects passing filter in question
0049     trigger::TriggerObjectCollection allTriggerObjects = triggerEvent->getObjects();  // all objects
0050     trigger::TriggerObjectCollection selectedObjects;  // vector to fill with objects passing our filter
0051     if (filterIndex < triggerEvent->sizeFilters()) {
0052       const auto& keys = triggerEvent->filterKeys(filterIndex);
0053       for (unsigned short key : keys) {
0054         trigger::TriggerObject foundObject = allTriggerObjects[key];
0055         selectedObjects.push_back(foundObject);
0056       }
0057     }
0058 
0059     // differentiate between event level and particle level variables
0060     const static std::vector<std::string> eventLevelVariables = {"AK4HT", "AK8HT", "MET"};
0061     if (std::find(eventLevelVariables.begin(), eventLevelVariables.end(), objType_) != eventLevelVariables.end()) {
0062       // for these event level variables we only require the existence of a trigger object, but no matching
0063       if (!selectedObjects.empty())
0064         for (auto& hist : hists_)
0065           hist->fill(obj);
0066     } else {
0067       // do a deltaR matching between trigger object and GEN object
0068       double mindR2 = 99999;
0069       for (const auto& filterobj : selectedObjects) {
0070         double dR = deltaR2(obj, filterobj);
0071         if (dR < mindR2)
0072           mindR2 = dR;
0073       }
0074       if (mindR2 < dR2limit_)
0075         for (auto& hist : hists_)
0076           hist->fill(obj);
0077     }
0078   }
0079 }
0080 
0081 // booker function for 1D hists
0082 void HLTGenValHistCollFilter::book1D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig) {
0083   // extracting parameters from configuration
0084   auto vsVar = histConfig.getParameter<std::string>("vsVar");
0085   auto vsVarFunc = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVar);
0086 
0087   // this class parses any potential additional cuts, changes in binning or tags
0088   HLTGenValPathSpecificSettingParser parser =
0089       HLTGenValPathSpecificSettingParser(histConfig.getParameter<std::string>("pathSpecificCuts"),
0090                                          histConfig.getParameter<std::vector<edm::ParameterSet>>("binnings"),
0091                                          vsVar);
0092 
0093   // getting the bin edges, path-specific overwrites general if present
0094   auto binLowEdgesDouble = histConfig.getParameter<std::vector<double>>("binLowEdges");
0095   if (parser.havePathSpecificBins())
0096     binLowEdgesDouble = *parser.getPathSpecificBins();
0097 
0098   // additional cuts applied to this histogram, combination of general ones and path-specific ones
0099   std::vector<edm::ParameterSet> allCutsVector = histConfig.getParameter<std::vector<edm::ParameterSet>>("rangeCuts");
0100   std::vector<edm::ParameterSet> pathSpecificCutsVector = *parser.getPathSpecificCuts();
0101   allCutsVector.insert(allCutsVector.end(), pathSpecificCutsVector.begin(), pathSpecificCutsVector.end());
0102   VarRangeCutColl<HLTGenValObject> rangeCuts(allCutsVector);
0103 
0104   // getting the custom tag
0105   const std::string tag = *parser.getTag();
0106 
0107   // checking validity of vsVar
0108   if (!vsVarFunc) {
0109     throw cms::Exception("ConfigError") << " vsVar " << vsVar << " is giving null ptr (likely empty) in " << __FILE__
0110                                         << "," << __LINE__ << std::endl;
0111   }
0112 
0113   // converting bin edges to float
0114   std::vector<float> binLowEdges;
0115   binLowEdges.reserve(binLowEdgesDouble.size());
0116   for (double lowEdge : binLowEdgesDouble)
0117     binLowEdges.push_back(lowEdge);
0118 
0119   // name and title are systematically constructed
0120 
0121   // remove potential leading "-" (which denotes that that trigger is ignored)
0122   std::string filterName = filter_;
0123   if (filterName.rfind('-', 0) == 0)
0124     filterName.erase(0, 1);
0125 
0126   std::string histName, histTitle;
0127   if (filter_ == "beforeAnyFilter") {  // this handles the naming of the "before" hist
0128     histName = objType_ + separator_ + path_ + separator_ + "GEN" + separator_ + "vs" + vsVar;
0129     histTitle = objType_ + " " + path_ + " GEN vs " + vsVar;
0130     if (!tag.empty()) {
0131       histName += separator_ + tag;
0132       histTitle += " " + tag;
0133     }
0134   } else {  // naming of all regular hists
0135     histName = objType_ + separator_ + path_ + separator_ + filterName + separator_ + "vs" + vsVar;
0136     histTitle = objType_ + " " + path_ + " " + filterName + " vs" + vsVar;
0137 
0138     // appending the tag, in case it is filled
0139     if (!tag.empty()) {
0140       histName += separator_ + tag;
0141       histTitle += " " + tag;
0142     }
0143   }
0144 
0145   auto me = iBooker.book1D(
0146       histName.c_str(), histTitle.c_str(), binLowEdges.size() - 1, binLowEdges.data());  // booking MonitorElement
0147 
0148   std::unique_ptr<HLTGenValHist> hist;  // creating the hist object
0149 
0150   hist = std::make_unique<HLTGenValHist1D>(me->getTH1(), vsVar, vsVarFunc, rangeCuts);
0151 
0152   hists_.emplace_back(std::move(hist));
0153 }
0154 
0155 // booker function for 2D hists
0156 void HLTGenValHistCollFilter::book2D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig2D) {
0157   // extracting parameters from configuration
0158   auto vsVarX = histConfig2D.getParameter<std::string>("vsVarX");
0159   auto vsVarY = histConfig2D.getParameter<std::string>("vsVarY");
0160   auto vsVarFuncX = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVarX);
0161   auto vsVarFuncY = hltdqm::getUnaryFuncFloat<HLTGenValObject>(vsVarY);
0162   auto binLowEdgesDoubleX = histConfig2D.getParameter<std::vector<double>>("binLowEdgesX");
0163   auto binLowEdgesDoubleY = histConfig2D.getParameter<std::vector<double>>("binLowEdgesY");
0164 
0165   // checking validity of vsVar
0166   if (!vsVarFuncX) {
0167     throw cms::Exception("ConfigError") << " vsVar " << vsVarX << " is giving null ptr (likely empty) in " << __FILE__
0168                                         << "," << __LINE__ << std::endl;
0169   }
0170   if (!vsVarFuncY) {
0171     throw cms::Exception("ConfigError") << " vsVar " << vsVarY << " is giving null ptr (likely empty) in " << __FILE__
0172                                         << "," << __LINE__ << std::endl;
0173   }
0174 
0175   // converting bin edges to float
0176   std::vector<float> binLowEdgesX;
0177   std::vector<float> binLowEdgesY;
0178   binLowEdgesX.reserve(binLowEdgesDoubleX.size());
0179   binLowEdgesY.reserve(binLowEdgesDoubleY.size());
0180   for (double lowEdge : binLowEdgesDoubleX)
0181     binLowEdgesX.push_back(lowEdge);
0182   for (double lowEdge : binLowEdgesDoubleY)
0183     binLowEdgesY.push_back(lowEdge);
0184 
0185   // name and title are systematically constructed
0186 
0187   // remove potential leading "-" (which denotes that that trigger is ignored)
0188   std::string filterName = filter_;
0189   if (filterName.rfind('-', 0) == 0)
0190     filterName.erase(0, 1);
0191 
0192   std::string histName, histTitle;
0193   if (filter_ == "beforeAnyFilter") {
0194     histName = objType_ + separator_ + path_ + separator_ + "GEN" + separator_ + "2Dvs" + vsVarX + separator_ + vsVarY;
0195     histTitle = objType_ + " " + path_ + " GEN 2D vs " + vsVarX + " " + vsVarY;
0196   } else {
0197     histName = objType_ + separator_ + path_ + separator_ + filterName + separator_ + "2Dvs" + vsVarX + vsVarY;
0198     histTitle = objType_ + " " + path_ + " " + filterName + " 2D vs" + vsVarX + " " + vsVarY;
0199   }
0200 
0201   auto me = iBooker.book2D(histName.c_str(),
0202                            histTitle.c_str(),
0203                            binLowEdgesX.size() - 1,
0204                            binLowEdgesX.data(),
0205                            binLowEdgesY.size() - 1,
0206                            binLowEdgesY.data());
0207 
0208   std::unique_ptr<HLTGenValHist> hist;
0209 
0210   hist = std::make_unique<HLTGenValHist2D>(me->getTH2F(), vsVarX, vsVarY, vsVarFuncX, vsVarFuncY);
0211 
0212   hists_.emplace_back(std::move(hist));
0213 }