Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:55

0001 /*
0002  *  Class:DQMGenericClient 
0003  *
0004  *  DQM histogram post processor
0005  *
0006  * 
0007  *  \author Junghwan Goh - SungKyunKwan University
0008  */
0009 
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ServiceRegistry/interface/Service.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0015 
0016 #include <TH1.h>
0017 #include <TH1F.h>
0018 #include <TH2F.h>
0019 #include <TClass.h>
0020 #include <TString.h>
0021 #include <TPRegexp.h>
0022 #include <TDirectory.h>
0023 #include <TEfficiency.h>
0024 
0025 #include <set>
0026 #include <cmath>
0027 #include <string>
0028 #include <vector>
0029 #include <climits>
0030 #include <boost/tokenizer.hpp>
0031 
0032 using namespace std;
0033 using namespace edm;
0034 
0035 class DQMGenericClient : public DQMEDHarvester {
0036 public:
0037   DQMGenericClient(const edm::ParameterSet& pset);
0038   ~DQMGenericClient() override {}
0039 
0040   void dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0041                              DQMStore::IGetter& igetter,
0042                              const edm::LuminosityBlock& lumiSeg,
0043                              const edm::EventSetup& c) override;
0044   void dqmEndRun(DQMStore::IBooker&, DQMStore::IGetter&, edm::Run const&, edm::EventSetup const&) override;
0045   void dqmEndJob(DQMStore::IBooker&, DQMStore::IGetter&) override {}
0046 
0047   enum class EfficType { none = 0, efficiency, fakerate, simpleratio };
0048 
0049   struct EfficOption {
0050     std::string name, title;
0051     std::string numerator, denominator;
0052     EfficType type;
0053     bool isProfile;
0054   };
0055 
0056   struct ResolOption {
0057     std::string namePrefix, titlePrefix;
0058     std::string srcName;
0059   };
0060 
0061   struct ProfileOption {
0062     std::string name, title;
0063     std::string srcName;
0064   };
0065 
0066   struct NormOption {
0067     std::string name, normHistName;
0068   };
0069 
0070   struct CDOption {
0071     std::string name;
0072     bool ascending;
0073   };
0074 
0075   struct NoFlowOption {
0076     std::string name;
0077   };
0078 
0079   void computeEfficiency(DQMStore::IBooker& ibooker,
0080                          DQMStore::IGetter& igetter,
0081                          const std::string& startDir,
0082                          const std::string& efficMEName,
0083                          const std::string& efficMETitle,
0084                          const std::string& recoMEName,
0085                          const std::string& simMEName,
0086                          const EfficType type = EfficType::efficiency,
0087                          const bool makeProfile = false);
0088   void computeResolution(DQMStore::IBooker& ibooker,
0089                          DQMStore::IGetter& igetter,
0090                          const std::string& startDir,
0091                          const std::string& fitMEPrefix,
0092                          const std::string& fitMETitlePrefix,
0093                          const std::string& srcMEName);
0094   void computeProfile(DQMStore::IBooker& ibooker,
0095                       DQMStore::IGetter& igetter,
0096                       const std::string& startDir,
0097                       const std::string& profileMEName,
0098                       const std::string& profileMETitle,
0099                       const std::string& srcMEName);
0100 
0101   void normalizeToEntries(DQMStore::IBooker& ibooker,
0102                           DQMStore::IGetter& igetter,
0103                           const std::string& startDir,
0104                           const std::string& histName,
0105                           const std::string& normHistName);
0106   void makeCumulativeDist(DQMStore::IBooker& ibooker,
0107                           DQMStore::IGetter& igetter,
0108                           const std::string& startDir,
0109                           const std::string& cdName,
0110                           bool ascending = true);
0111   void makeNoFlowDist(DQMStore::IBooker& ibooker,
0112                       DQMStore::IGetter& igetter,
0113                       const std::string& startDir,
0114                       const std::string& cdName);
0115 
0116   void limitedFit(MonitorElement* srcME, MonitorElement* meanME, MonitorElement* sigmaME);
0117 
0118 private:
0119   TPRegexp metacharacters_;
0120   TPRegexp nonPerlWildcard_;
0121   unsigned int verbose_;
0122   bool runOnEndLumi_;
0123   bool runOnEndJob_;
0124   bool makeGlobalEffPlot_;
0125   bool isWildcardUsed_;
0126   bool resLimitedFit_;
0127 
0128   DQMStore* theDQM;
0129   std::vector<std::string> subDirs_;
0130   std::string outputFileName_;
0131 
0132   std::vector<EfficOption> efficOptions_;
0133   std::vector<ResolOption> resolOptions_;
0134   std::vector<ProfileOption> profileOptions_;
0135   std::vector<NormOption> normOptions_;
0136   std::vector<CDOption> cdOptions_;
0137   std::vector<NoFlowOption> noFlowOptions_;
0138 
0139   void generic_eff(TH1* denom, TH1* numer, MonitorElement* efficiencyHist, const EfficType type = EfficType::efficiency);
0140 
0141   void findAllSubdirectories(DQMStore::IBooker& ibooker,
0142                              DQMStore::IGetter& igetter,
0143                              std::string dir,
0144                              std::set<std::string>* myList,
0145                              const TString& pattern);
0146 
0147   void makeAllPlots(DQMStore::IBooker&, DQMStore::IGetter&);
0148 
0149   void removeMEIfBooked(const std::string& meName, DQMStore::IGetter& igetter);
0150 };
0151 
0152 class FitSlicesYTool {
0153 public:
0154   typedef dqm::harvesting::MonitorElement MonitorElement;
0155   FitSlicesYTool(MonitorElement* me) {
0156     const bool oldAddDir = TH1::AddDirectoryStatus();
0157     TH1::AddDirectory(true);
0158     // ... create your hists
0159     TH2F* h = me->getTH2F();
0160     TF1 fgaus("fgaus", "gaus", h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax(), TF1::EAddToList::kNo);
0161     h->FitSlicesY(&fgaus, 0, -1, 0, "QNR SERIAL");
0162     string name(h->GetName());
0163     h0 = (TH1*)gDirectory->Get((name + "_0").c_str());
0164     h1 = (TH1*)gDirectory->Get((name + "_1").c_str());
0165     h2 = (TH1*)gDirectory->Get((name + "_2").c_str());
0166     h3 = (TH1*)gDirectory->Get((name + "_chi2").c_str());
0167     TH1::AddDirectory(oldAddDir);
0168   }
0169 
0170   /// Destructor
0171   ~FitSlicesYTool() {
0172     delete h0;
0173     delete h1;
0174     delete h2;
0175     delete h3;
0176   }
0177   /// Fill the ME with the mean value of the gaussian fit in each slice
0178   void getFittedMean(MonitorElement* me) {
0179     if (!(h1 && me))
0180       throw cms::Exception("FitSlicesYTool") << "Pointer =0 : h1=" << h1 << " me=" << me;
0181     if (h1->GetNbinsX() == me->getNbinsX()) {
0182       for (int bin = 0; bin != h1->GetNbinsX(); bin++) {
0183         me->setBinContent(bin + 1, h1->GetBinContent(bin + 1));
0184         //       me->setBinEntries(bin+1, 1.);
0185       }
0186     } else {
0187       throw cms::Exception("FitSlicesYTool") << "Different number of bins!";
0188     }
0189   }
0190   /// Fill the ME with the sigma value of the gaussian fit in each slice
0191   void getFittedSigma(MonitorElement* me) {
0192     if (!(h2 && me))
0193       throw cms::Exception("FitSlicesYTool") << "Pointer =0 : h1=" << h1 << " me=" << me;
0194     if (h2->GetNbinsX() == me->getNbinsX()) {
0195       for (int bin = 0; bin != h2->GetNbinsX(); bin++) {
0196         me->setBinContent(bin + 1, h2->GetBinContent(bin + 1));
0197         //       me->setBinEntries(bin+1, 1.);
0198       }
0199     } else {
0200       throw cms::Exception("FitSlicesYTool") << "Different number of bins!";
0201     }
0202   }
0203   /// Fill the ME with the mean value (with error) of the gaussian fit in each slice
0204   void getFittedMeanWithError(MonitorElement* me) {
0205     if (!(h1 && me))
0206       throw cms::Exception("FitSlicesYTool") << "Pointer =0 : h1=" << h1 << " me=" << me;
0207     if (h1->GetNbinsX() == me->getNbinsX()) {
0208       for (int bin = 0; bin != h1->GetNbinsX(); bin++) {
0209         me->setBinContent(bin + 1, h1->GetBinContent(bin + 1));
0210         //       me->setBinEntries(bin+1, 1.);
0211         me->setBinError(bin + 1, h1->GetBinError(bin + 1));
0212       }
0213     } else {
0214       throw cms::Exception("FitSlicesYTool") << "Different number of bins!";
0215     }
0216   }
0217   /// Fill the ME with the sigma value (with error) of the gaussian fit in each slice
0218   void getFittedSigmaWithError(MonitorElement* me) {
0219     if (!(h2 && me))
0220       throw cms::Exception("FitSlicesYTool") << "Pointer =0 : h1=" << h1 << " me=" << me;
0221     if (h2->GetNbinsX() == me->getNbinsX()) {
0222       for (int bin = 0; bin != h2->GetNbinsX(); bin++) {
0223         me->setBinContent(bin + 1, h2->GetBinContent(bin + 1));
0224         //       me->setBinEntries(bin+1, 1.);
0225         me->setBinError(bin + 1, h2->GetBinError(bin + 1));
0226       }
0227     } else {
0228       throw cms::Exception("FitSlicesYTool") << "Different number of bins!";
0229     }
0230   }
0231 
0232 private:
0233   TH1* h0;
0234   TH1* h1;
0235   TH1* h2;
0236   TH1* h3;
0237 };
0238 
0239 typedef DQMGenericClient::MonitorElement ME;
0240 
0241 DQMGenericClient::DQMGenericClient(const ParameterSet& pset)
0242     : metacharacters_("[\\^\\$\\.\\*\\+\\?\\|\\(\\)\\{\\}\\[\\]]"), nonPerlWildcard_("\\w\\*|^\\*") {
0243   typedef std::vector<edm::ParameterSet> VPSet;
0244   typedef std::vector<std::string> vstring;
0245   typedef boost::escaped_list_separator<char> elsc;
0246 
0247   elsc commonEscapes("\\", " \t", "\'");
0248 
0249   verbose_ = pset.getUntrackedParameter<unsigned int>("verbose", 0);
0250   runOnEndLumi_ = pset.getUntrackedParameter<bool>("runOnEndLumi", false);
0251   runOnEndJob_ = pset.getUntrackedParameter<bool>("runOnEndJob", true);
0252   makeGlobalEffPlot_ = pset.getUntrackedParameter<bool>("makeGlobalEffienciesPlot", true);
0253 
0254   // Parse efficiency commands
0255   vstring effCmds = pset.getParameter<vstring>("efficiency");
0256   for (vstring::const_iterator effCmd = effCmds.begin(); effCmd != effCmds.end(); ++effCmd) {
0257     if (effCmd->empty())
0258       continue;
0259 
0260     boost::tokenizer<elsc> tokens(*effCmd, commonEscapes);
0261 
0262     vector<string> args;
0263     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0264       if (iToken->empty())
0265         continue;
0266       args.push_back(*iToken);
0267     }
0268 
0269     if (args.size() < 4) {
0270       LogInfo("DQMGenericClient") << "Wrong input to effCmds\n";
0271       continue;
0272     }
0273 
0274     EfficOption opt;
0275     opt.name = args[0];
0276     opt.title = args[1];
0277     opt.numerator = args[2];
0278     opt.denominator = args[3];
0279     opt.isProfile = false;
0280 
0281     const string typeName = args.size() == 4 ? "eff" : args[4];
0282     if (typeName == "eff")
0283       opt.type = EfficType::efficiency;
0284     else if (typeName == "fake")
0285       opt.type = EfficType::fakerate;
0286     else if (typeName == "simpleratio")
0287       opt.type = EfficType::simpleratio;
0288     else
0289       opt.type = EfficType::none;
0290 
0291     efficOptions_.push_back(opt);
0292   }
0293 
0294   VPSet efficSets = pset.getUntrackedParameter<VPSet>("efficiencySets", VPSet());
0295   for (VPSet::const_iterator efficSet = efficSets.begin(); efficSet != efficSets.end(); ++efficSet) {
0296     EfficOption opt;
0297     opt.name = efficSet->getUntrackedParameter<string>("name");
0298     opt.title = efficSet->getUntrackedParameter<string>("title");
0299     opt.numerator = efficSet->getUntrackedParameter<string>("numerator");
0300     opt.denominator = efficSet->getUntrackedParameter<string>("denominator");
0301     opt.isProfile = false;
0302 
0303     const string typeName = efficSet->getUntrackedParameter<string>("typeName", "eff");
0304     if (typeName == "eff")
0305       opt.type = EfficType::efficiency;
0306     else if (typeName == "fake")
0307       opt.type = EfficType::fakerate;
0308     else if (typeName == "simpleratio")
0309       opt.type = EfficType::simpleratio;
0310     else
0311       opt.type = EfficType::none;
0312 
0313     efficOptions_.push_back(opt);
0314   }
0315 
0316   // Parse efficiency profiles
0317   vstring effProfileCmds = pset.getUntrackedParameter<vstring>("efficiencyProfile", vstring());
0318   for (vstring::const_iterator effProfileCmd = effProfileCmds.begin(); effProfileCmd != effProfileCmds.end();
0319        ++effProfileCmd) {
0320     if (effProfileCmd->empty())
0321       continue;
0322 
0323     boost::tokenizer<elsc> tokens(*effProfileCmd, commonEscapes);
0324 
0325     vector<string> args;
0326     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0327       if (iToken->empty())
0328         continue;
0329       args.push_back(*iToken);
0330     }
0331 
0332     if (args.size() < 4) {
0333       LogInfo("DQMGenericClient") << "Wrong input to effProfileCmds\n";
0334       continue;
0335     }
0336 
0337     EfficOption opt;
0338     opt.name = args[0];
0339     opt.title = args[1];
0340     opt.numerator = args[2];
0341     opt.denominator = args[3];
0342     opt.isProfile = true;
0343 
0344     const string typeName = args.size() == 4 ? "eff" : args[4];
0345     if (typeName == "eff")
0346       opt.type = EfficType::efficiency;
0347     else if (typeName == "fake")
0348       opt.type = EfficType::fakerate;
0349     else if (typeName == "simpleratio")
0350       opt.type = EfficType::simpleratio;
0351     else
0352       opt.type = EfficType::none;
0353 
0354     efficOptions_.push_back(opt);
0355   }
0356 
0357   VPSet effProfileSets = pset.getUntrackedParameter<VPSet>("efficiencyProfileSets", VPSet());
0358   for (VPSet::const_iterator effProfileSet = effProfileSets.begin(); effProfileSet != effProfileSets.end();
0359        ++effProfileSet) {
0360     EfficOption opt;
0361     opt.name = effProfileSet->getUntrackedParameter<string>("name");
0362     opt.title = effProfileSet->getUntrackedParameter<string>("title");
0363     opt.numerator = effProfileSet->getUntrackedParameter<string>("numerator");
0364     opt.denominator = effProfileSet->getUntrackedParameter<string>("denominator");
0365     opt.isProfile = true;
0366 
0367     const string typeName = effProfileSet->getUntrackedParameter<string>("typeName", "eff");
0368     if (typeName == "eff")
0369       opt.type = EfficType::efficiency;
0370     else if (typeName == "fake")
0371       opt.type = EfficType::fakerate;
0372     else if (typeName == "simpleratio")
0373       opt.type = EfficType::simpleratio;
0374     else
0375       opt.type = EfficType::none;
0376 
0377     efficOptions_.push_back(opt);
0378   }
0379 
0380   // Parse resolution commands
0381   vstring resCmds = pset.getParameter<vstring>("resolution");
0382   for (vstring::const_iterator resCmd = resCmds.begin(); resCmd != resCmds.end(); ++resCmd) {
0383     if (resCmd->empty())
0384       continue;
0385     boost::tokenizer<elsc> tokens(*resCmd, commonEscapes);
0386 
0387     vector<string> args;
0388     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0389       if (iToken->empty())
0390         continue;
0391       args.push_back(*iToken);
0392     }
0393 
0394     if (args.size() != 3) {
0395       LogInfo("DQMGenericClient") << "Wrong input to resCmds\n";
0396       continue;
0397     }
0398 
0399     ResolOption opt;
0400     opt.namePrefix = args[0];
0401     opt.titlePrefix = args[1];
0402     opt.srcName = args[2];
0403 
0404     resolOptions_.push_back(opt);
0405   }
0406 
0407   VPSet resolSets = pset.getUntrackedParameter<VPSet>("resolutionSets", VPSet());
0408   for (VPSet::const_iterator resolSet = resolSets.begin(); resolSet != resolSets.end(); ++resolSet) {
0409     ResolOption opt;
0410     opt.namePrefix = resolSet->getUntrackedParameter<string>("namePrefix");
0411     opt.titlePrefix = resolSet->getUntrackedParameter<string>("titlePrefix");
0412     opt.srcName = resolSet->getUntrackedParameter<string>("srcName");
0413 
0414     resolOptions_.push_back(opt);
0415   }
0416 
0417   // Parse profiles
0418   vstring profileCmds = pset.getUntrackedParameter<vstring>("profile", vstring());
0419   for (const auto& profileCmd : profileCmds) {
0420     boost::tokenizer<elsc> tokens(profileCmd, commonEscapes);
0421 
0422     vector<string> args;
0423     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0424       if (iToken->empty())
0425         continue;
0426       args.push_back(*iToken);
0427     }
0428 
0429     if (args.size() != 3) {
0430       LogInfo("DQMGenericClient") << "Wrong input to profileCmds\n";
0431       continue;
0432     }
0433 
0434     ProfileOption opt;
0435     opt.name = args[0];
0436     opt.title = args[1];
0437     opt.srcName = args[2];
0438 
0439     profileOptions_.push_back(opt);
0440   }
0441 
0442   VPSet profileSets = pset.getUntrackedParameter<VPSet>("profileSets", VPSet());
0443   for (const auto& profileSet : profileSets) {
0444     ProfileOption opt;
0445     opt.name = profileSet.getUntrackedParameter<string>("name");
0446     opt.title = profileSet.getUntrackedParameter<string>("title");
0447     opt.srcName = profileSet.getUntrackedParameter<string>("srcName");
0448 
0449     profileOptions_.push_back(opt);
0450   }
0451 
0452   // Parse Normalization commands
0453   vstring normCmds = pset.getUntrackedParameter<vstring>("normalization", vstring());
0454   for (vstring::const_iterator normCmd = normCmds.begin(); normCmd != normCmds.end(); ++normCmd) {
0455     if (normCmd->empty())
0456       continue;
0457     boost::tokenizer<elsc> tokens(*normCmd, commonEscapes);
0458 
0459     vector<string> args;
0460     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0461       if (iToken->empty())
0462         continue;
0463       args.push_back(*iToken);
0464     }
0465 
0466     if (args.empty() or args.size() > 2) {
0467       LogInfo("DQMGenericClient") << "Wrong input to normCmds\n";
0468       continue;
0469     }
0470 
0471     NormOption opt;
0472     opt.name = args[0];
0473     opt.normHistName = args.size() == 2 ? args[1] : args[0];
0474 
0475     normOptions_.push_back(opt);
0476   }
0477 
0478   VPSet normSets = pset.getUntrackedParameter<VPSet>("normalizationSets", VPSet());
0479   for (VPSet::const_iterator normSet = normSets.begin(); normSet != normSets.end(); ++normSet) {
0480     NormOption opt;
0481     opt.name = normSet->getUntrackedParameter<string>("name");
0482     opt.normHistName = normSet->getUntrackedParameter<string>("normalizedTo", opt.name);
0483 
0484     normOptions_.push_back(opt);
0485   }
0486 
0487   // Cumulative distributions
0488   vstring cdCmds = pset.getUntrackedParameter<vstring>("cumulativeDists", vstring());
0489   for (vstring::const_iterator cdCmd = cdCmds.begin(); cdCmd != cdCmds.end(); ++cdCmd) {
0490     if (cdCmd->empty())
0491       continue;
0492     boost::tokenizer<elsc> tokens(*cdCmd, commonEscapes);
0493 
0494     vector<string> args;
0495     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0496       if (iToken->empty())
0497         continue;
0498       args.push_back(*iToken);
0499     }
0500 
0501     if (args.empty() || args.size() > 2) {
0502       LogInfo("DQMGenericClient") << "Wrong input to cdCmds\n";
0503       continue;
0504     }
0505 
0506     CDOption opt;
0507     opt.name = args[0];
0508     opt.ascending = args.size() == 2 ? (args[1] != "descending") : true;
0509 
0510     cdOptions_.push_back(opt);
0511   }
0512 
0513   VPSet cdSets = pset.getUntrackedParameter<VPSet>("cumulativeDistSets", VPSet());
0514   for (VPSet::const_iterator cdSet = cdSets.begin(); cdSet != cdSets.end(); ++cdSet) {
0515     CDOption opt;
0516     opt.name = cdSet->getUntrackedParameter<string>("name");
0517     opt.ascending = cdSet->getUntrackedParameter<bool>("ascending", true);
0518 
0519     cdOptions_.push_back(opt);
0520   }
0521 
0522   // move under/overflows to first/last bins
0523   vstring noFlowCmds = pset.getUntrackedParameter<vstring>("noFlowDists", vstring());
0524   for (vstring::const_iterator noFlowCmd = noFlowCmds.begin(); noFlowCmd != noFlowCmds.end(); ++noFlowCmd) {
0525     if (noFlowCmd->empty())
0526       continue;
0527     boost::tokenizer<elsc> tokens(*noFlowCmd, commonEscapes);
0528 
0529     vector<string> args;
0530     for (boost::tokenizer<elsc>::const_iterator iToken = tokens.begin(); iToken != tokens.end(); ++iToken) {
0531       if (iToken->empty())
0532         continue;
0533       args.push_back(*iToken);
0534     }
0535 
0536     if (args.empty() || args.size() > 2) {
0537       LogInfo("DQMGenericClient") << "Wrong input to noFlowCmds\n";
0538       continue;
0539     }
0540 
0541     NoFlowOption opt;
0542     opt.name = args[0];
0543 
0544     noFlowOptions_.push_back(opt);
0545   }
0546 
0547   VPSet noFlowSets = pset.getUntrackedParameter<VPSet>("noFlowDistSets", VPSet());
0548   for (VPSet::const_iterator noFlowSet = noFlowSets.begin(); noFlowSet != noFlowSets.end(); ++noFlowSet) {
0549     NoFlowOption opt;
0550     opt.name = noFlowSet->getUntrackedParameter<string>("name");
0551 
0552     noFlowOptions_.push_back(opt);
0553   }
0554 
0555   outputFileName_ = pset.getUntrackedParameter<string>("outputFileName", "");
0556   subDirs_ = pset.getUntrackedParameter<vstring>("subDirs");
0557 
0558   resLimitedFit_ = pset.getUntrackedParameter<bool>("resolutionLimitedFit", false);
0559   isWildcardUsed_ = false;
0560 }
0561 
0562 void DQMGenericClient::dqmEndLuminosityBlock(DQMStore::IBooker& ibooker,
0563                                              DQMStore::IGetter& igetter,
0564                                              const edm::LuminosityBlock& lumiSeg,
0565                                              const edm::EventSetup& c) {
0566   if (runOnEndLumi_) {
0567     makeAllPlots(ibooker, igetter);
0568   }
0569 }
0570 
0571 void DQMGenericClient::dqmEndRun(DQMStore::IBooker& ibooker,
0572                                  DQMStore::IGetter& igetter,
0573                                  edm::Run const&,
0574                                  edm::EventSetup const&) {
0575   // Create new MEs in endRun, even though we are requested to do it in endJob.
0576   // This gives the QTests a chance to run, before summaries are created in
0577   // endJob. The negative side effect is that we cannot run the GenericClient
0578   // for plots produced in Harvesting, but that seems rather rare.
0579   //
0580   // It is important that this is still save in the presence of multiple runs,
0581   // first because in multi-run harvesting, we accumulate statistics over all
0582   // runs and have full statistics at the endRun of the last run, and second,
0583   // because we set the efficiencyFlag so any further aggregation should produce
0584   // correct results. Also, all operations should be idempotent; running them
0585   // more than once does no harm.
0586 
0587   // needed to access the DQMStore::save method
0588   theDQM = nullptr;
0589   theDQM = Service<DQMStore>().operator->();
0590 
0591   if (runOnEndJob_) {
0592     makeAllPlots(ibooker, igetter);
0593   }
0594 
0595   if (!outputFileName_.empty())
0596     theDQM->save(outputFileName_);
0597 }
0598 
0599 void DQMGenericClient::makeAllPlots(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0600   typedef vector<string> vstring;
0601 
0602   // Process wildcard in the sub-directory
0603   set<string> subDirSet;
0604 
0605   for (vstring::const_iterator iSubDir = subDirs_.begin(); iSubDir != subDirs_.end(); ++iSubDir) {
0606     string subDir = *iSubDir;
0607 
0608     if (subDir[subDir.size() - 1] == '/')
0609       subDir.erase(subDir.size() - 1);
0610 
0611     if (TString(subDir).Contains(metacharacters_)) {
0612       isWildcardUsed_ = true;
0613 
0614       const string::size_type shiftPos = subDir.rfind('/');
0615       const string searchPath = subDir.substr(0, shiftPos);
0616       const string pattern = subDir.substr(shiftPos + 1, subDir.length());
0617       //std::cout << "\n\n\n\nLooking for all subdirs of " << subDir << std::endl;
0618 
0619       findAllSubdirectories(ibooker, igetter, searchPath, &subDirSet, pattern);
0620 
0621     } else {
0622       subDirSet.insert(subDir);
0623     }
0624   }
0625 
0626   for (set<string>::const_iterator iSubDir = subDirSet.begin(); iSubDir != subDirSet.end(); ++iSubDir) {
0627     const string& dirName = *iSubDir;
0628 
0629     // First normalize, then move under/overflows, then make
0630     // cumulative, and only then efficiency This allows to use the
0631     // cumulative distributions for efficiency calculation
0632     for (vector<NormOption>::const_iterator normOption = normOptions_.begin(); normOption != normOptions_.end();
0633          ++normOption) {
0634       normalizeToEntries(ibooker, igetter, dirName, normOption->name, normOption->normHistName);
0635     }
0636 
0637     for (vector<NoFlowOption>::const_iterator noFlowOption = noFlowOptions_.begin();
0638          noFlowOption != noFlowOptions_.end();
0639          ++noFlowOption) {
0640       makeNoFlowDist(ibooker, igetter, dirName, noFlowOption->name);
0641     }
0642 
0643     for (vector<CDOption>::const_iterator cdOption = cdOptions_.begin(); cdOption != cdOptions_.end(); ++cdOption) {
0644       makeCumulativeDist(ibooker, igetter, dirName, cdOption->name, cdOption->ascending);
0645     }
0646 
0647     for (vector<EfficOption>::const_iterator efficOption = efficOptions_.begin(); efficOption != efficOptions_.end();
0648          ++efficOption) {
0649       computeEfficiency(ibooker,
0650                         igetter,
0651                         dirName,
0652                         efficOption->name,
0653                         efficOption->title,
0654                         efficOption->numerator,
0655                         efficOption->denominator,
0656                         efficOption->type,
0657                         efficOption->isProfile);
0658     }
0659     if (makeGlobalEffPlot_) {
0660       ME* globalEfficME = igetter.get(dirName + "/globalEfficiencies");
0661       if (globalEfficME)
0662         globalEfficME->getTH1F()->LabelsDeflate("X");
0663     }
0664 
0665     for (vector<ResolOption>::const_iterator resolOption = resolOptions_.begin(); resolOption != resolOptions_.end();
0666          ++resolOption) {
0667       computeResolution(
0668           ibooker, igetter, dirName, resolOption->namePrefix, resolOption->titlePrefix, resolOption->srcName);
0669     }
0670 
0671     for (const auto& profileOption : profileOptions_) {
0672       computeProfile(ibooker, igetter, dirName, profileOption.name, profileOption.title, profileOption.srcName);
0673     }
0674   }
0675 }
0676 
0677 void DQMGenericClient::computeEfficiency(DQMStore::IBooker& ibooker,
0678                                          DQMStore::IGetter& igetter,
0679                                          const string& startDir,
0680                                          const string& efficMEName,
0681                                          const string& efficMETitle,
0682                                          const string& recoMEName,
0683                                          const string& simMEName,
0684                                          const EfficType type,
0685                                          const bool makeProfile) {
0686   if (!igetter.dirExists(startDir)) {
0687     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0688       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0689                                   << "Cannot find sub-directory " << startDir << endl;
0690     }
0691     return;
0692   }
0693 
0694   ibooker.cd();
0695 
0696   ME* simME = igetter.get(startDir + "/" + simMEName);
0697   ME* recoME = igetter.get(startDir + "/" + recoMEName);
0698 
0699   if (!simME) {
0700     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0701       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0702                                   << "No sim-ME '" << simMEName << "' found\n";
0703     }
0704     return;
0705   }
0706 
0707   if (!recoME) {
0708     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0709       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0710                                   << "No reco-ME '" << recoMEName << "' found\n";
0711     }
0712     return;
0713   }
0714 
0715   // Treat everything as the base class, TH1
0716 
0717   TH1* hSim = simME->getTH1();
0718   TH1* hReco = recoME->getTH1();
0719 
0720   if (!hSim || !hReco) {
0721     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0722       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0723                                   << "Cannot create TH1 from ME\n";
0724     }
0725     return;
0726   }
0727 
0728   string efficDir = startDir;
0729   string newEfficMEName = efficMEName;
0730   string::size_type shiftPos;
0731   if (string::npos != (shiftPos = efficMEName.rfind('/'))) {
0732     efficDir += "/" + efficMEName.substr(0, shiftPos);
0733     newEfficMEName.erase(0, shiftPos + 1);
0734   }
0735   ibooker.setCurrentFolder(efficDir);
0736 
0737   if (makeProfile) {
0738     TProfile* efficHist = (hReco->GetXaxis()->GetXbins()->GetSize() == 0)
0739                               ? new TProfile(newEfficMEName.c_str(),
0740                                              efficMETitle.c_str(),
0741                                              hReco->GetXaxis()->GetNbins(),
0742                                              hReco->GetXaxis()->GetXmin(),
0743                                              hReco->GetXaxis()->GetXmax())
0744                               : new TProfile(newEfficMEName.c_str(),
0745                                              efficMETitle.c_str(),
0746                                              hReco->GetXaxis()->GetNbins(),
0747                                              hReco->GetXaxis()->GetXbins()->GetArray());
0748 
0749     efficHist->GetXaxis()->SetTitle(hSim->GetXaxis()->GetTitle());
0750     efficHist->GetYaxis()->SetTitle(hSim->GetYaxis()->GetTitle());
0751 
0752     for (int i = 1; i <= hReco->GetNbinsX(); i++) {
0753       const double nReco = hReco->GetBinContent(i);
0754       const double nSim = hSim->GetBinContent(i);
0755 
0756       if (!std::string(hSim->GetXaxis()->GetBinLabel(i)).empty())
0757         efficHist->GetXaxis()->SetBinLabel(i, hSim->GetXaxis()->GetBinLabel(i));
0758 
0759       if (nSim == 0 or nReco < 0 or nReco > nSim)
0760         continue;
0761       const double effVal = nReco / nSim;
0762       const double errLo = TEfficiency::ClopperPearson(nSim, nReco, 0.683, false);
0763       const double errUp = TEfficiency::ClopperPearson(nSim, nReco, 0.683, true);
0764       const double errVal = (effVal - errLo > errUp - effVal) ? effVal - errLo : errUp - effVal;
0765       efficHist->SetBinContent(i, effVal);
0766       efficHist->SetBinEntries(i, 1);
0767       efficHist->SetBinError(i, std::hypot(effVal, errVal));
0768     }
0769     ibooker.bookProfile(newEfficMEName, efficHist);
0770     delete efficHist;
0771   }
0772 
0773   else {
0774     TH1* efficHist = (TH1*)hSim->Clone(newEfficMEName.c_str());
0775     efficHist->SetTitle(efficMETitle.c_str());
0776 
0777     // Here is where you have trouble --- you need
0778     // to understand what type of hist you have.
0779 
0780     ME* efficME = nullptr;
0781 
0782     // Parse the class name
0783     // This works, but there might be a better way
0784     TClass* myHistClass = efficHist->IsA();
0785     TString histClassName = myHistClass->GetName();
0786 
0787     if (histClassName == "TH1F") {
0788       efficME = ibooker.book1D(newEfficMEName, (TH1F*)efficHist);
0789     } else if (histClassName == "TH2F") {
0790       efficME = ibooker.book2D(newEfficMEName, (TH2F*)efficHist);
0791     } else if (histClassName == "TH3F") {
0792       efficME = ibooker.book3D(newEfficMEName, (TH3F*)efficHist);
0793     }
0794 
0795     delete efficHist;
0796 
0797     if (!efficME) {
0798       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0799                                   << "Cannot book effic-ME from the DQM\n";
0800       return;
0801     }
0802 
0803     // Update: 2009-9-16 slaunwhj
0804     // call the most generic efficiency function
0805     // works up to 3-d histograms
0806 
0807     generic_eff(hSim, hReco, efficME, type);
0808 
0809     //   const int nBin = efficME->getNbinsX();
0810     //   for(int bin = 0; bin <= nBin; ++bin) {
0811     //     const float nSim  = simME ->getBinContent(bin);
0812     //     const float nReco = recoME->getBinContent(bin);
0813     //     float eff =0;
0814     //     if (type=="fake")eff = nSim ? 1-nReco/nSim : 0.;
0815     //     else eff= nSim ? nReco/nSim : 0.;
0816     //     const float err = nSim && eff <= 1 ? sqrt(eff*(1-eff)/nSim) : 0.;
0817     //     efficME->setBinContent(bin, eff);
0818     //     efficME->setBinError(bin, err);
0819     //   }
0820     efficME->setEntries(simME->getEntries());
0821   }
0822 
0823   // Global efficiency
0824   if (makeGlobalEffPlot_) {
0825     ME* globalEfficME = igetter.get(efficDir + "/globalEfficiencies");
0826     if (!globalEfficME)
0827       globalEfficME = ibooker.book1D("globalEfficiencies", "Global efficiencies", 1, 0, 1);
0828     if (!globalEfficME) {
0829       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0830                                   << "Cannot book globalEffic-ME from the DQM\n";
0831       return;
0832     }
0833     globalEfficME->setEfficiencyFlag();
0834     TH1F* hGlobalEffic = globalEfficME->getTH1F();
0835     if (!hGlobalEffic) {
0836       LogInfo("DQMGenericClient") << "computeEfficiency() : "
0837                                   << "Cannot create TH1F from ME, globalEfficME\n";
0838       return;
0839     }
0840 
0841     const float nSimAll = hSim->GetEntries();
0842     const float nRecoAll = hReco->GetEntries();
0843     float efficAll = 0;
0844     if (type == EfficType::efficiency || type == EfficType::simpleratio)
0845       efficAll = nSimAll ? nRecoAll / nSimAll : 0;
0846     else if (type == EfficType::fakerate)
0847       efficAll = nSimAll ? 1 - nRecoAll / nSimAll : 0;
0848     float errorAll = 0;
0849     if (type == EfficType::simpleratio) {
0850       if (nSimAll) {
0851         const float x = nRecoAll / nSimAll;
0852         errorAll = std::sqrt(1.f / nSimAll * x * (1 + x));
0853       }
0854     } else
0855       errorAll = nSimAll && efficAll < 1 ? sqrt(efficAll * (1 - efficAll) / nSimAll) : 0;
0856 
0857     const int iBin = hGlobalEffic->Fill(newEfficMEName.c_str(), 0);
0858     hGlobalEffic->SetBinContent(iBin, efficAll);
0859     hGlobalEffic->SetBinError(iBin, errorAll);
0860   }
0861 }
0862 
0863 void DQMGenericClient::computeResolution(DQMStore::IBooker& ibooker,
0864                                          DQMStore::IGetter& igetter,
0865                                          const string& startDir,
0866                                          const string& namePrefix,
0867                                          const string& titlePrefix,
0868                                          const std::string& srcName) {
0869   if (!igetter.dirExists(startDir)) {
0870     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0871       LogInfo("DQMGenericClient") << "computeResolution() : "
0872                                   << "Cannot find sub-directory " << startDir << endl;
0873     }
0874     return;
0875   }
0876 
0877   ibooker.cd();
0878 
0879   ME* srcME = igetter.get(startDir + "/" + srcName);
0880   if (!srcME) {
0881     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0882       LogInfo("DQMGenericClient") << "computeResolution() : "
0883                                   << "No source ME '" << srcName << "' found\n";
0884     }
0885     return;
0886   }
0887 
0888   TH2F* hSrc = srcME->getTH2F();
0889   if (!hSrc) {
0890     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0891       LogInfo("DQMGenericClient") << "computeResolution() : "
0892                                   << "Cannot create TH2F from source-ME\n";
0893     }
0894     return;
0895   }
0896 
0897   const int nBin = hSrc->GetNbinsX();
0898 
0899   string newDir = startDir;
0900   string newPrefix = namePrefix;
0901   string::size_type shiftPos;
0902   if (string::npos != (shiftPos = namePrefix.rfind('/'))) {
0903     newDir += "/" + namePrefix.substr(0, shiftPos);
0904     newPrefix.erase(0, shiftPos + 1);
0905   }
0906 
0907   ibooker.setCurrentFolder(newDir);
0908 
0909   float* lowedgesfloats = new float[nBin + 1];
0910   ME* meanME;
0911   ME* sigmaME;
0912   if (hSrc->GetXaxis()->GetXbins()->GetSize()) {
0913     for (int j = 0; j < nBin + 1; ++j)
0914       lowedgesfloats[j] = (float)hSrc->GetXaxis()->GetXbins()->GetAt(j);
0915     meanME = ibooker.book1D(newPrefix + "_Mean", titlePrefix + " Mean", nBin, lowedgesfloats);
0916     sigmaME = ibooker.book1D(newPrefix + "_Sigma", titlePrefix + " Sigma", nBin, lowedgesfloats);
0917   } else {
0918     meanME = ibooker.book1D(
0919         newPrefix + "_Mean", titlePrefix + " Mean", nBin, hSrc->GetXaxis()->GetXmin(), hSrc->GetXaxis()->GetXmax());
0920     sigmaME = ibooker.book1D(
0921         newPrefix + "_Sigma", titlePrefix + " Sigma", nBin, hSrc->GetXaxis()->GetXmin(), hSrc->GetXaxis()->GetXmax());
0922   }
0923 
0924   if (meanME && sigmaME) {
0925     meanME->setEfficiencyFlag();
0926     sigmaME->setEfficiencyFlag();
0927 
0928     if (!resLimitedFit_) {
0929       FitSlicesYTool fitTool(srcME);
0930       fitTool.getFittedMeanWithError(meanME);
0931       fitTool.getFittedSigmaWithError(sigmaME);
0932       ////  fitTool.getFittedChisqWithError(chi2ME); // N/A
0933     } else {
0934       limitedFit(srcME, meanME, sigmaME);
0935     }
0936   }
0937   delete[] lowedgesfloats;
0938 }
0939 
0940 void DQMGenericClient::computeProfile(DQMStore::IBooker& ibooker,
0941                                       DQMStore::IGetter& igetter,
0942                                       const std::string& startDir,
0943                                       const std::string& profileMEName,
0944                                       const std::string& profileMETitle,
0945                                       const std::string& srcMEName) {
0946   if (!igetter.dirExists(startDir)) {
0947     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0948       LogInfo("DQMGenericClient") << "computeProfile() : "
0949                                   << "Cannot find sub-directory " << startDir << endl;
0950     }
0951     return;
0952   }
0953 
0954   ibooker.cd();
0955 
0956   ME* srcME = igetter.get(startDir + "/" + srcMEName);
0957   if (!srcME) {
0958     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0959       LogInfo("DQMGenericClient") << "computeProfile() : "
0960                                   << "No source ME '" << srcMEName << "' found\n";
0961     }
0962     return;
0963   }
0964 
0965   TH2F* hSrc = srcME->getTH2F();
0966   if (!hSrc) {
0967     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0968       LogInfo("DQMGenericClient") << "computeProfile() : "
0969                                   << "Cannot create TH2F from source-ME\n";
0970     }
0971     return;
0972   }
0973 
0974   string profileDir = startDir;
0975   string newProfileMEName = profileMEName;
0976   string::size_type shiftPos;
0977   if (string::npos != (shiftPos = profileMEName.rfind('/'))) {
0978     profileDir += "/" + profileMEName.substr(0, shiftPos);
0979     newProfileMEName.erase(0, shiftPos + 1);
0980   }
0981   ibooker.setCurrentFolder(profileDir);
0982 
0983   std::unique_ptr<TProfile> profile(hSrc->ProfileX());  // We own the pointer
0984   profile->SetTitle(profileMETitle.c_str());
0985   ibooker.bookProfile(profileMEName, profile.get());  // ibooker makes a copy
0986 }
0987 
0988 void DQMGenericClient::normalizeToEntries(DQMStore::IBooker& ibooker,
0989                                           DQMStore::IGetter& igetter,
0990                                           const std::string& startDir,
0991                                           const std::string& histName,
0992                                           const std::string& normHistName) {
0993   if (!igetter.dirExists(startDir)) {
0994     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
0995       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
0996                                   << "Cannot find sub-directory " << startDir << endl;
0997     }
0998     return;
0999   }
1000 
1001   ibooker.cd();
1002 
1003   ME* element = igetter.get(startDir + "/" + histName);
1004   ME* normME = igetter.get(startDir + "/" + normHistName);
1005 
1006   if (!element) {
1007     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1008       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
1009                                   << "No such element '" << histName << "' found\n";
1010     }
1011     return;
1012   }
1013 
1014   if (!normME) {
1015     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1016       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
1017                                   << "No such element '" << normHistName << "' found\n";
1018     }
1019     return;
1020   }
1021 
1022   TH1F* hist = element->getTH1F();
1023   if (!hist) {
1024     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1025       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
1026                                   << "Cannot create TH1F from ME\n";
1027     }
1028     return;
1029   }
1030 
1031   TH1F* normHist = normME->getTH1F();
1032   if (!normHist) {
1033     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1034       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
1035                                   << "Cannot create TH1F from ME\n";
1036     }
1037     return;
1038   }
1039 
1040   const double entries = normHist->GetEntries();
1041   if (entries != 0) {
1042     hist->Scale(1. / entries);
1043   } else {
1044     LogInfo("DQMGenericClient") << "normalizeToEntries() : "
1045                                 << "Zero entries in histogram\n";
1046   }
1047 
1048   return;
1049 }
1050 
1051 void DQMGenericClient::makeCumulativeDist(DQMStore::IBooker& ibooker,
1052                                           DQMStore::IGetter& igetter,
1053                                           const std::string& startDir,
1054                                           const std::string& cdName,
1055                                           bool ascending) {
1056   if (!igetter.dirExists(startDir)) {
1057     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1058       LogInfo("DQMGenericClient") << "makeCumulativeDist() : "
1059                                   << "Cannot find sub-directory " << startDir << endl;
1060     }
1061     return;
1062   }
1063 
1064   ibooker.cd();
1065 
1066   ME* element_cd = igetter.get(startDir + "/" + cdName);
1067 
1068   if (!element_cd) {
1069     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1070       LogInfo("DQMGenericClient") << "makeCumulativeDist() : "
1071                                   << "No such element '" << cdName << "' found\n";
1072     }
1073     return;
1074   }
1075 
1076   TH1F* cd = element_cd->getTH1F();
1077 
1078   if (!cd) {
1079     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1080       LogInfo("DQMGenericClient") << "makeCumulativeDist() : "
1081                                   << "Cannot create TH1F from ME\n";
1082     }
1083     return;
1084   }
1085 
1086   int n_bins = cd->GetNbinsX() + 1;
1087 
1088   if (ascending) {
1089     for (int i = 1; i <= n_bins; i++) {
1090       cd->SetBinContent(i, cd->GetBinContent(i) + cd->GetBinContent(i - 1));
1091     }
1092   } else {
1093     for (int i = n_bins - 1; i >= 0; i--) {  // n_bins points to the overflow bin
1094       cd->SetBinContent(i, cd->GetBinContent(i) + cd->GetBinContent(i + 1));
1095     }
1096   }
1097 
1098   return;
1099 }
1100 
1101 void DQMGenericClient::makeNoFlowDist(DQMStore::IBooker& ibooker,
1102                                       DQMStore::IGetter& igetter,
1103                                       const std::string& startDir,
1104                                       const std::string& noFlowName) {
1105   if (!igetter.dirExists(startDir)) {
1106     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1107       LogInfo("DQMGenericClient") << "makeNoFlowDist() : "
1108                                   << "Cannot find sub-directory " << startDir << endl;
1109     }
1110     return;
1111   }
1112 
1113   ibooker.cd();
1114 
1115   ME* element_noFlow = igetter.get(startDir + "/" + noFlowName);
1116 
1117   if (!element_noFlow) {
1118     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1119       LogInfo("DQMGenericClient") << "makeNoFlowDist() : "
1120                                   << "No such element '" << noFlowName << "' found\n";
1121     }
1122     return;
1123   }
1124 
1125   TH1F* noFlow = element_noFlow->getTH1F();
1126 
1127   if (!noFlow) {
1128     if (verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_)) {
1129       LogInfo("DQMGenericClient") << "makeNoFlowDist() : "
1130                                   << "Cannot create TH1F from ME\n";
1131     }
1132     return;
1133   }
1134 
1135   noFlow->AddBinContent(1, noFlow->GetBinContent(0));
1136   noFlow->SetBinContent(0, 0.);
1137 
1138   const auto lastBin = noFlow->GetNbinsX();
1139   noFlow->AddBinContent(lastBin, noFlow->GetBinContent(lastBin + 1));
1140   noFlow->SetBinContent(lastBin + 1, 0.);
1141 }
1142 
1143 void DQMGenericClient::limitedFit(MonitorElement* srcME, MonitorElement* meanME, MonitorElement* sigmaME) {
1144   TH2F* histo = srcME->getTH2F();
1145 
1146   // Fit slices projected along Y from bins in X
1147   double cont_min = 100;  //Minimum number of entries
1148   Int_t binx = histo->GetXaxis()->GetNbins();
1149 
1150   for (int i = 1; i <= binx; i++) {
1151     TString iString(i);
1152     TH1* histoY = histo->ProjectionY(" ", i, i);
1153     double cont = histoY->GetEntries();
1154 
1155     if (cont >= cont_min) {
1156       float minfit = histoY->GetMean() - histoY->GetRMS();
1157       float maxfit = histoY->GetMean() + histoY->GetRMS();
1158 
1159       TF1* fitFcn = new TF1(TString("g") + histo->GetName() + iString, "gaus", minfit, maxfit);
1160       double x1, x2;
1161       fitFcn->GetRange(x1, x2);
1162       //To avoid DQMGenericClient maintains state between fit calls
1163       fitFcn->SetParameters(histoY->Integral(), histoY->GetMean(), histoY->GetRMS());
1164 
1165       histoY->Fit(fitFcn, "QR0 SERIAL", "", x1, x2);
1166 
1167       //      histoY->Fit(fitFcn->GetName(),"RME");
1168       double* par = fitFcn->GetParameters();
1169       const double* err = fitFcn->GetParErrors();
1170 
1171       meanME->setBinContent(i, par[1]);
1172       meanME->setBinError(i, err[1]);
1173       //       meanME->setBinEntries(i, 1.);
1174       //       meanME->setBinError(i,sqrt(err[1]*err[1]+par[1]*par[1]));
1175 
1176       sigmaME->setBinContent(i, par[2]);
1177       sigmaME->setBinError(i, err[2]);
1178       //       sigmaME->setBinEntries(i, 1.);
1179       //       sigmaME->setBinError(i,sqrt(err[2]*err[2]+par[2]*par[2]));
1180 
1181       if (fitFcn)
1182         delete fitFcn;
1183       if (histoY)
1184         delete histoY;
1185     } else {
1186       if (histoY)
1187         delete histoY;
1188       continue;
1189     }
1190   }
1191 }
1192 
1193 //=================================
1194 
1195 void DQMGenericClient::findAllSubdirectories(DQMStore::IBooker& ibooker,
1196                                              DQMStore::IGetter& igetter,
1197                                              std::string dir,
1198                                              std::set<std::string>* myList,
1199                                              const TString& _pattern = TString("")) {
1200   TString pattern = _pattern;
1201   if (!igetter.dirExists(dir)) {
1202     LogError("DQMGenericClient") << " DQMGenericClient::findAllSubdirectories ==> Missing folder " << dir << " !!!";
1203     return;
1204   }
1205   if (pattern != "") {
1206     if (pattern.Contains(nonPerlWildcard_))
1207       pattern.ReplaceAll("*", ".*");
1208     TPRegexp regexp(pattern);
1209     ibooker.cd(dir);
1210     vector<string> foundDirs = igetter.getSubdirs();
1211     for (vector<string>::const_iterator iDir = foundDirs.begin(); iDir != foundDirs.end(); ++iDir) {
1212       TString dirName = iDir->substr(iDir->rfind('/') + 1, iDir->length());
1213       if (dirName.Contains(regexp))
1214         findAllSubdirectories(ibooker, igetter, *iDir, myList);
1215     }
1216   }
1217   //std::cout << "Looking for directory " << dir ;
1218   else if (igetter.dirExists(dir)) {
1219     //std::cout << "... it exists! Inserting it into the list ";
1220     myList->insert(dir);
1221     //std::cout << "... now list has size " << myList->size() << std::endl;
1222     ibooker.cd(dir);
1223     findAllSubdirectories(ibooker, igetter, dir, myList, "*");
1224   } else {
1225     //std::cout << "... DOES NOT EXIST!!! Skip bogus dir" << std::endl;
1226 
1227     LogInfo("DQMGenericClient") << "Trying to find sub-directories of " << dir << " failed because " << dir
1228                                 << " does not exist";
1229   }
1230   return;
1231 }
1232 
1233 void DQMGenericClient::generic_eff(TH1* denom, TH1* numer, MonitorElement* efficiencyHist, const EfficType type) {
1234   for (int iBinX = 1; iBinX < denom->GetNbinsX() + 1; iBinX++) {
1235     for (int iBinY = 1; iBinY < denom->GetNbinsY() + 1; iBinY++) {
1236       for (int iBinZ = 1; iBinZ < denom->GetNbinsZ() + 1; iBinZ++) {
1237         int globalBinNum = denom->GetBin(iBinX, iBinY, iBinZ);
1238 
1239         float numerVal = numer->GetBinContent(globalBinNum);
1240         float denomVal = denom->GetBinContent(globalBinNum);
1241 
1242         float effVal = 0;
1243 
1244         // fake eff is in use
1245         if (type == EfficType::fakerate) {
1246           effVal = denomVal ? (1 - numerVal / denomVal) : 0;
1247         } else {
1248           effVal = denomVal ? numerVal / denomVal : 0;
1249         }
1250 
1251         float errVal = 0;
1252         if (type == EfficType::simpleratio) {
1253           //          errVal = denomVal ? 1.f/denomVal*effVal*(1+effVal) : 0;
1254           float numerErr = numer->GetBinError(globalBinNum);
1255           float denomErr = denom->GetBinError(globalBinNum);
1256           float denomsq = denomVal * denomVal;
1257           errVal = denomVal ? sqrt(pow(1.f / denomVal * numerErr, 2.0) + pow(numerVal / denomsq * denomErr, 2)) : 0;
1258         } else {
1259           errVal = (denomVal && (effVal <= 1)) ? sqrt(effVal * (1 - effVal) / denomVal) : 0;
1260         }
1261 
1262         LogDebug("DQMGenericClient") << "(iBinX, iBinY, iBinZ)  = " << iBinX << ", " << iBinY << ", " << iBinZ
1263                                      << "), global bin =  " << globalBinNum << "eff = " << numerVal << "  /  "
1264                                      << denomVal << " =  " << effVal << " ... setting the error for that bin ... "
1265                                      << endl
1266                                      << endl;
1267 
1268         efficiencyHist->setBinContent(globalBinNum, effVal);
1269         efficiencyHist->setBinError(globalBinNum, errVal);
1270         efficiencyHist->setEfficiencyFlag();
1271       }
1272     }
1273   }
1274 
1275   //efficiencyHist->setMinimum(0.0);
1276   //efficiencyHist->setMaximum(1.0);
1277 }
1278 
1279 DEFINE_FWK_MODULE(DQMGenericClient);
1280 
1281 /* vim:set ts=2 sts=2 sw=2 expandtab: */