DQMGenericTnPClient

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DQMServices/Core/interface/DQMStore.h"

#include "DQMOffline/Trigger/plugins/GenericTnPFitter.h"

#include <TString.h>
#include <TPRegexp.h>

using namespace edm;
using namespace dqmTnP;
using namespace std;

using vstring = std::vector<std::string>;

class DQMGenericTnPClient : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
public:
  typedef dqm::legacy::MonitorElement MonitorElement;
  typedef dqm::legacy::DQMStore DQMStore;

  DQMGenericTnPClient(const edm::ParameterSet& pset);
  ~DQMGenericTnPClient() override;
  void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override {}
  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
  void endRun(const edm::Run& run, const edm::EventSetup& setup) override;
  void calculateEfficiency(const std::string& dirName, const ParameterSet& pset);
  void findAllSubdirectories(const std::string& dir, std::set<std::string>* myList, TString pattern);

private:
  DQMStore* dqmStore;
  TFile* plots;
  vstring subDirs;
  std::string myDQMrootFolder;
  bool verbose;
  const VParameterSet efficiencies;
  GaussianPlusLinearFitter* GPLfitter;
  VoigtianPlusExponentialFitter* VPEfitter;
};

DQMGenericTnPClient::DQMGenericTnPClient(const edm::ParameterSet& pset)
    : subDirs(pset.getUntrackedParameter<vstring>("subDirs", vstring())),
      myDQMrootFolder(pset.getUntrackedParameter<std::string>("MyDQMrootFolder", "")),
      verbose(pset.getUntrackedParameter<bool>("Verbose", false)),
      efficiencies(pset.getUntrackedParameter<VParameterSet>("Efficiencies")) {
  usesResource("DQMStore");
  TString savePlotsInRootFileName = pset.getUntrackedParameter<string>("SavePlotsInRootFileName", "");
  plots = savePlotsInRootFileName != "" ? new TFile(savePlotsInRootFileName, "recreate") : nullptr;
  GPLfitter = new GaussianPlusLinearFitter(verbose);
  VPEfitter = new VoigtianPlusExponentialFitter(verbose);
}

void DQMGenericTnPClient::endRun(const edm::Run& run, const edm::EventSetup& setup) {
  TPRegexp metacharacters("[\\^\\$\\.\\*\\+\\?\\|\\(\\)\\{\\}\\[\\]]");

  dqmStore = Service<DQMStore>().operator->();
  if (!dqmStore) {
    LogError("DQMGenericTnPClient") << "Could not find DQMStore service\n";
    return;
  }
  dqmStore->setCurrentFolder(myDQMrootFolder);

  set<std::string> subDirSet;

  if (!myDQMrootFolder.empty())
    subDirSet.insert(myDQMrootFolder);
  else {
    for (auto iSubDir = subDirs.begin(); iSubDir != subDirs.end(); ++iSubDir) {
      std::string subDir = *iSubDir;
      if (subDir[subDir.size() - 1] == '/')
        subDir.erase(subDir.size() - 1);
      if (TString(subDir).Contains(metacharacters)) {
        const string::size_type shiftPos = subDir.rfind('/');
        const string searchPath = subDir.substr(0, shiftPos);
        const string pattern = subDir.substr(shiftPos + 1, subDir.length());
        findAllSubdirectories(searchPath, &subDirSet, pattern);
      } else {
        subDirSet.insert(subDir);
      }
    }
  }

  for (auto const& dirName : subDirSet) {
    for (auto const& efficiencie : efficiencies) {
      calculateEfficiency(dirName, efficiencie);
    }
  }
}

void DQMGenericTnPClient::calculateEfficiency(const std::string& dirName, const ParameterSet& pset) {
  //get hold of numerator and denominator histograms
  string allMEname = dirName + "/" + pset.getUntrackedParameter<string>("DenominatorMEname");
  string passMEname = dirName + "/" + pset.getUntrackedParameter<string>("NumeratorMEname");
  MonitorElement* allME = dqmStore->get(allMEname);
  MonitorElement* passME = dqmStore->get(passMEname);
  if (allME == nullptr || passME == nullptr) {
    LogDebug("DQMGenericTnPClient") << "Could not find MEs: " << allMEname << " or " << passMEname << endl;
    return;
  }
  TH1* all = allME->getTH1();
  TH1* pass = passME->getTH1();
  //setup the fitter
  string fitFunction = pset.getUntrackedParameter<string>("FitFunction");
  AbstractFitter* fitter = nullptr;
  if (fitFunction == "GaussianPlusLinear") {
    GPLfitter->setup(pset.getUntrackedParameter<double>("ExpectedMean"),
                     pset.getUntrackedParameter<double>("FitRangeLow"),
                     pset.getUntrackedParameter<double>("FitRangeHigh"),
                     pset.getUntrackedParameter<double>("ExpectedSigma"));
    fitter = GPLfitter;
  } else if (fitFunction == "VoigtianPlusExponential") {
    VPEfitter->setup(pset.getUntrackedParameter<double>("ExpectedMean"),
                     pset.getUntrackedParameter<double>("FitRangeLow"),
                     pset.getUntrackedParameter<double>("FitRangeHigh"),
                     pset.getUntrackedParameter<double>("ExpectedSigma"),
                     pset.getUntrackedParameter<double>("FixedWidth"));
    fitter = VPEfitter;
  } else {
    LogError("DQMGenericTnPClient") << "Fit function: " << fitFunction << " does not exist" << endl;
    return;
  }
  //check dimensions
  int dimensions = all->GetDimension();
  int massDimension = pset.getUntrackedParameter<int>("MassDimension");
  if (massDimension > dimensions) {
    LogError("DQMGenericTnPClient") << "Monitoring Element " << allMEname << " has smaller dimension than "
                                    << massDimension << endl;
    return;
  }
  //figure out directory and efficiency names
  string effName = pset.getUntrackedParameter<string>("EfficiencyMEname");
  string effDir = dirName;
  string::size_type slashPos = effName.rfind('/');
  if (string::npos != slashPos) {
    effDir += "/" + effName.substr(0, slashPos);
    effName.erase(0, slashPos + 1);
  }
  dqmStore->setCurrentFolder(effDir);
  TString prefix(effDir.c_str());
  prefix.ReplaceAll('/', '_');
  //calculate and book efficiency
  if (dimensions == 2) {
    TProfile* eff = nullptr;
    TProfile* effChi2 = nullptr;
    TString error = fitter->calculateEfficiency(
        (TH2*)pass, (TH2*)all, massDimension, eff, effChi2, plots ? prefix + effName.c_str() : "");
    if (error != "") {
      LogError("DQMGenericTnPClient") << error << endl;
      return;
    }
    dqmStore->bookProfile(effName, eff);
    dqmStore->bookProfile(effName + "Chi2", effChi2);
    delete eff;
    delete effChi2;
  } else if (dimensions == 3) {
    TProfile2D* eff = nullptr;
    TProfile2D* effChi2 = nullptr;
    TString error = fitter->calculateEfficiency(
        (TH3*)pass, (TH3*)all, massDimension, eff, effChi2, plots ? prefix + effName.c_str() : "");
    if (error != "") {
      LogError("DQMGenericTnPClient") << error << endl;
      return;
    }
    dqmStore->bookProfile2D(effName, eff);
    dqmStore->bookProfile2D(effName + "Chi2", effChi2);
    delete eff;
    delete effChi2;
  }
}

DQMGenericTnPClient::~DQMGenericTnPClient() {
  delete GPLfitter;
  if (plots) {
    plots->Close();
  }
}

void DQMGenericTnPClient::findAllSubdirectories(const std::string& dir,
                                                std::set<std::string>* myList,
                                                TString pattern = "") {
  if (!dqmStore->dirExists(dir)) {
    LogError("DQMGenericTnPClient") << " DQMGenericTnPClient::findAllSubdirectories ==> Missing folder " << dir
                                    << " !!!";
    return;
  }
  TPRegexp nonPerlWildcard("\\w\\*|^\\*");
  if (pattern != "") {
    if (pattern.Contains(nonPerlWildcard))
      pattern.ReplaceAll("*", ".*");
    TPRegexp regexp(pattern);
    dqmStore->cd(dir);
    vector<string> foundDirs = dqmStore->getSubdirs();
    for (auto iDir = foundDirs.begin(); iDir != foundDirs.end(); ++iDir) {
      TString dirName = iDir->substr(iDir->rfind('/') + 1, iDir->length());
      if (dirName.Contains(regexp))
        findAllSubdirectories(*iDir, myList);
    }
  } else if (dqmStore->dirExists(dir)) {
    myList->insert(dir);
    dqmStore->cd(dir);
    findAllSubdirectories(dir, myList, "*");
  } else {
    LogInfo("DQMGenericClient") << "Trying to find sub-directories of " << dir << " failed because " << dir
                                << " does not exist";
  }
  return;
}

//define this as a plug-in
DEFINE_FWK_MODULE(DQMGenericTnPClient);