Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:12

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/ServiceRegistry/interface/Service.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 
0009 #include "PhysicsTools/TagAndProbe/interface/TagProbeFitter.h"
0010 
0011 using namespace std;
0012 using namespace edm;
0013 
0014 class TagProbeFitTreeAnalyzer : public edm::one::EDAnalyzer<> {
0015 public:
0016   TagProbeFitTreeAnalyzer(const edm::ParameterSet& pset);
0017   ~TagProbeFitTreeAnalyzer() override{};
0018   void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override{};
0019   void calculateEfficiency(string name, const edm::ParameterSet& pset);
0020 
0021 private:
0022   TagProbeFitter fitter;
0023   unsigned int split_mode;  // number of events to read per cycle (slower, but memory efficient)
0024 };
0025 
0026 TagProbeFitTreeAnalyzer::TagProbeFitTreeAnalyzer(const edm::ParameterSet& pset)
0027     : fitter(
0028           pset.getParameter<vector<string> >("InputFileNames"),
0029           pset.getParameter<string>("InputDirectoryName"),
0030           pset.getParameter<string>("InputTreeName"),
0031           pset.getParameter<string>("OutputFileName"),
0032           pset.existsAs<unsigned int>("NumCPU") ? pset.getParameter<unsigned int>("NumCPU") : 1,
0033           pset.existsAs<bool>("SaveWorkspace") ? pset.getParameter<bool>("SaveWorkspace") : false,
0034           pset.existsAs<bool>("floatShapeParameters") ? pset.getParameter<bool>("floatShapeParameters") : true,
0035           pset.existsAs<vector<string> >("fixVars") ? pset.getParameter<vector<string> >("fixVars") : vector<string>()),
0036       split_mode(pset.existsAs<unsigned int>("SplitMode") ? pset.getParameter<unsigned int>("SplitMode") : 0) {
0037   fitter.setQuiet(pset.getUntrackedParameter("Quiet", false));
0038   fitter.setSplitMode(split_mode);
0039 
0040   if (pset.existsAs<bool>("binnedFit")) {
0041     bool binned = pset.getParameter<bool>("binnedFit");
0042     fitter.setBinnedFit(binned, binned ? pset.getParameter<uint32_t>("binsForFit") : 0);
0043   } else if (pset.existsAs<uint32_t>("binsForMassPlots")) {
0044     fitter.setBinsForMassPlots(pset.getParameter<uint32_t>("binsForMassPlots"));
0045   }
0046 
0047   if (pset.existsAs<bool>("saveDistributionsPlot")) {
0048     fitter.setSaveDistributionsPlot(pset.getParameter<bool>("saveDistributionsPlot"));
0049   }
0050   if (pset.existsAs<std::string>("WeightVariable")) {
0051     fitter.setWeightVar(pset.getParameter<std::string>("WeightVariable"));
0052   }
0053   const ParameterSet variables = pset.getParameter<ParameterSet>("Variables");
0054   vector<string> variableNames = variables.getParameterNamesForType<vector<string> >();
0055   for (vector<string>::const_iterator name = variableNames.begin(); name != variableNames.end(); name++) {
0056     vector<string> var = variables.getParameter<vector<string> >(*name);
0057     double lo, hi;
0058     if (var.size() >= 4 && !(istringstream(var[1]) >> lo).fail() && !(istringstream(var[2]) >> hi).fail()) {
0059       fitter.addVariable(*name, var[0], lo, hi, var[3]);
0060     } else {
0061       LogError("TagProbeFitTreeAnalyzer")
0062           << "Could not create variable: " << *name
0063           << ". Example: pt = cms.vstring(\"Probe pT\", \"1.0\", \"100.0\", \"GeV/c\") ";
0064     }
0065   }
0066 
0067   const ParameterSet categories = pset.getParameter<ParameterSet>("Categories");
0068   vector<string> categoryNames = categories.getParameterNamesForType<vector<string> >();
0069   for (vector<string>::const_iterator name = categoryNames.begin(); name != categoryNames.end(); name++) {
0070     vector<string> cat = categories.getParameter<vector<string> >(*name);
0071     if (cat.size() == 2) {
0072       fitter.addCategory(*name, cat[0], cat[1]);
0073     } else {
0074       LogError("TagProbeFitTreeAnalyzer") << "Could not create category: " << *name
0075                                           << ". Example: mcTrue = cms.vstring(\"MC True\", \"dummy[true=1,false=0]\") ";
0076     }
0077   }
0078 
0079   if (pset.existsAs<ParameterSet>("Expressions")) {
0080     const ParameterSet exprs = pset.getParameter<ParameterSet>("Expressions");
0081     vector<string> exprNames = exprs.getParameterNamesForType<vector<string> >();
0082     for (vector<string>::const_iterator name = exprNames.begin(); name != exprNames.end(); name++) {
0083       vector<string> expr = exprs.getParameter<vector<string> >(*name);
0084       if (expr.size() >= 2) {
0085         vector<string> args(expr.begin() + 2, expr.end());
0086         fitter.addExpression(*name, expr[0], expr[1], args);
0087       } else {
0088         LogError("TagProbeFitTreeAnalyzer")
0089             << "Could not create expr: " << *name
0090             << ". Example: qop = cms.vstring(\"qOverP\", \"charge/p\", \"charge\", \"p\") ";
0091       }
0092     }
0093   }
0094 
0095   if (pset.existsAs<ParameterSet>("Cuts")) {
0096     const ParameterSet cuts = pset.getParameter<ParameterSet>("Cuts");
0097     vector<string> cutNames = cuts.getParameterNamesForType<vector<string> >();
0098     for (vector<string>::const_iterator name = cutNames.begin(); name != cutNames.end(); name++) {
0099       vector<string> cat = cuts.getParameter<vector<string> >(*name);
0100       if (cat.size() == 3) {
0101         fitter.addThresholdCategory(*name, cat[0], cat[1], atof(cat[2].c_str()));
0102       } else {
0103         LogError("TagProbeFitTreeAnalyzer") << "Could not create cut: " << *name
0104                                             << ". Example: matched = cms.vstring(\"Matched\", \"deltaR\", \"0.5\") ";
0105       }
0106     }
0107   }
0108 
0109   if (pset.existsAs<ParameterSet>("PDFs")) {
0110     const ParameterSet pdfs = pset.getParameter<ParameterSet>("PDFs");
0111     vector<string> pdfNames = pdfs.getParameterNamesForType<vector<string> >();
0112     for (vector<string>::const_iterator name = pdfNames.begin(); name != pdfNames.end(); name++) {
0113       vector<string> pdf = pdfs.getParameter<vector<string> >(*name);
0114       fitter.addPdf(*name, pdf);
0115     }
0116   }
0117 
0118   const ParameterSet efficiencies = pset.getParameter<ParameterSet>("Efficiencies");
0119   vector<string> efficiencyNames = efficiencies.getParameterNamesForType<ParameterSet>();
0120   for (vector<string>::const_iterator name = efficiencyNames.begin(); name != efficiencyNames.end(); name++) {
0121     try {
0122       calculateEfficiency(*name, efficiencies.getParameter<ParameterSet>(*name));
0123     } catch (std::exception& ex) {
0124       throw cms::Exception("Error", ex.what());
0125     }
0126   }
0127 }
0128 
0129 void TagProbeFitTreeAnalyzer::calculateEfficiency(string name, const edm::ParameterSet& pset) {
0130   vector<string> effCatState = pset.getParameter<vector<string> >("EfficiencyCategoryAndState");
0131   if (effCatState.empty() || (effCatState.size() % 2 == 1)) {
0132     cout << "EfficiencyCategoryAndState must be a even-sized list of category names and states of that category (cat1, "
0133             "state1, cat2, state2, ...)."
0134          << endl;
0135     exit(1);
0136   }
0137 
0138   vector<string> unbinnedVariables;
0139   if (pset.existsAs<vector<string> >("UnbinnedVariables")) {
0140     unbinnedVariables = pset.getParameter<vector<string> >("UnbinnedVariables");
0141   }
0142 
0143   const ParameterSet binVars = pset.getParameter<ParameterSet>("BinnedVariables");
0144   map<string, vector<double> > binnedVariables;
0145   vector<string> variableNames = binVars.getParameterNamesForType<vector<double> >();
0146   for (vector<string>::const_iterator var = variableNames.begin(); var != variableNames.end(); var++) {
0147     vector<double> binning = binVars.getParameter<vector<double> >(*var);
0148     binnedVariables[*var] = binning;
0149   }
0150   map<string, vector<string> > mappedCategories;
0151   vector<string> categoryNames = binVars.getParameterNamesForType<vector<string> >();
0152   for (vector<string>::const_iterator var = categoryNames.begin(); var != categoryNames.end(); var++) {
0153     vector<string> map = binVars.getParameter<vector<string> >(*var);
0154     mappedCategories[*var] = map;
0155   }
0156 
0157   vector<string> binToPDFmap;
0158   if (pset.existsAs<vector<string> >("BinToPDFmap")) {
0159     binToPDFmap = pset.getParameter<vector<string> >("BinToPDFmap");
0160   }
0161   if ((!binToPDFmap.empty()) && (binToPDFmap.size() % 2 == 0)) {
0162     cout << "BinToPDFmap must have odd size, first string is the default, followed by binRegExp - PDFname pairs!"
0163          << endl;
0164     exit(2);
0165   }
0166 
0167   vector<string> effCats, effStates;
0168   for (size_t i = 0, n = effCatState.size() / 2; i < n; ++i) {
0169     effCats.push_back(effCatState[2 * i]);
0170     effStates.push_back(effCatState[2 * i + 1]);
0171   }
0172 
0173   fitter.calculateEfficiency(
0174       name, effCats, effStates, unbinnedVariables, binnedVariables, mappedCategories, binToPDFmap);
0175 }
0176 
0177 //define this as a plug-in
0178 DEFINE_FWK_MODULE(TagProbeFitTreeAnalyzer);